Warning: this documentation for the development version is under construction.
00001 // Copyright (C) 2010 Marc Duruflé 00002 // 00003 // This file is part of the linear-algebra library Seldon, 00004 // http://seldon.sourceforge.net/. 00005 // 00006 // Seldon is free software; you can redistribute it and/or modify it under the 00007 // terms of the GNU Lesser General Public License as published by the Free 00008 // Software Foundation; either version 2.1 of the License, or (at your option) 00009 // any later version. 00010 // 00011 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY 00012 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 00013 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for 00014 // more details. 00015 // 00016 // You should have received a copy of the GNU Lesser General Public License 00017 // along with Seldon. If not, see http://www.gnu.org/licenses/. 00018 00019 00020 #ifndef SELDON_FILE_ILUT_PRECONDITIONING_CXX 00021 00022 #include "SymmetricIlutPreconditioning.cxx" 00023 00024 namespace Seldon 00025 { 00026 00027 template<class real, class cplx, class Allocator> 00028 IlutPreconditioning<real, cplx, Allocator>::IlutPreconditioning() 00029 { 00030 print_level = 0; 00031 symmetric_algorithm = false; 00032 type_ilu = ILUT; 00033 fill_level = 1000000; 00034 additional_fill = 1000000; 00035 mbloc = 1000000; 00036 alpha = 1.0; 00037 droptol = 0.01; 00038 permtol = 0.1; 00039 } 00040 00041 00042 template<class real, class cplx, class Allocator> 00043 void IlutPreconditioning<real, cplx, Allocator>::Clear() 00044 { 00045 permutation_row.Clear(); 00046 mat_sym.Clear(); 00047 mat_unsym.Clear(); 00048 xtmp.Clear(); 00049 } 00050 00051 00052 template<class real, class cplx, class Allocator> 00053 int IlutPreconditioning<real, cplx, Allocator>::GetFactorisationType() const 00054 { 00055 return type_ilu; 00056 } 00057 00058 00059 template<class real, class cplx, class Allocator> 00060 int IlutPreconditioning<real, cplx, Allocator>::GetFillLevel() const 00061 { 00062 return fill_level; 00063 } 00064 00065 00066 template<class real, class cplx, class Allocator> 00067 int IlutPreconditioning<real, cplx, Allocator>::GetAdditionalFillNumber() 00068 const 00069 { 00070 return additional_fill; 00071 } 00072 00073 00074 template<class real, class cplx, class Allocator> 00075 int IlutPreconditioning<real, cplx, Allocator>::GetPrintLevel() const 00076 { 00077 return print_level; 00078 } 00079 00080 00081 template<class real, class cplx, class Allocator> 00082 int IlutPreconditioning<real, cplx, Allocator>::GetPivotBlockInteger() const 00083 { 00084 return mbloc; 00085 } 00086 00087 00088 template<class real, class cplx, class Allocator> 00089 void IlutPreconditioning<real, cplx, Allocator> 00090 ::SetFactorisationType(int type) 00091 { 00092 type_ilu = type; 00093 } 00094 00095 00096 template<class real, class cplx, class Allocator> 00097 void IlutPreconditioning<real, cplx, Allocator>::SetFillLevel(int level) 00098 { 00099 fill_level = level; 00100 } 00101 00102 00103 template<class real, class cplx, class Allocator> 00104 void IlutPreconditioning<real, cplx, Allocator> 00105 ::SetAdditionalFillNumber(int level) 00106 { 00107 additional_fill = level; 00108 } 00109 00110 00111 template<class real, class cplx, class Allocator> 00112 void IlutPreconditioning<real, cplx, Allocator>::SetPrintLevel(int level) 00113 { 00114 print_level = level; 00115 } 00116 00117 00118 template<class real, class cplx, class Allocator> 00119 void IlutPreconditioning<real, cplx, Allocator>::SetPivotBlockInteger(int i) 00120 { 00121 mbloc = i; 00122 } 00123 00124 00125 template<class real, class cplx, class Allocator> 00126 real IlutPreconditioning<real, cplx, Allocator> 00127 ::GetDroppingThreshold() const 00128 { 00129 return droptol; 00130 } 00131 00132 00133 template<class real, class cplx, class Allocator> 00134 real IlutPreconditioning<real, cplx, Allocator> 00135 ::GetDiagonalCoefficient() const 00136 { 00137 return alpha; 00138 } 00139 00140 00141 template<class real, class cplx, class Allocator> 00142 real IlutPreconditioning<real, cplx, Allocator>::GetPivotThreshold() const 00143 { 00144 return permtol; 00145 } 00146 00147 00148 template<class real, class cplx, class Allocator> 00149 void IlutPreconditioning<real, cplx, Allocator> 00150 ::SetDroppingThreshold(real tol) 00151 { 00152 droptol = tol; 00153 } 00154 00155 00156 template<class real, class cplx, class Allocator> 00157 void IlutPreconditioning<real, cplx, Allocator> 00158 ::SetDiagonalCoefficient(real beta) 00159 { 00160 alpha = beta; 00161 } 00162 00163 00164 template<class real, class cplx, class Allocator> 00165 void IlutPreconditioning<real, cplx, Allocator>::SetPivotThreshold(real tol) 00166 { 00167 permtol = tol; 00168 } 00169 00170 00171 template<class real, class cplx, class Allocator> 00172 void IlutPreconditioning<real, cplx, Allocator>::SetSymmetricAlgorithm() 00173 { 00174 symmetric_algorithm = true; 00175 } 00176 00177 00178 template<class real, class cplx, class Allocator> 00179 void IlutPreconditioning<real, cplx, Allocator>::SetUnsymmetricAlgorithm() 00180 { 00181 symmetric_algorithm = false; 00182 } 00183 00184 00185 template<class real, class cplx, class Allocator> 00186 template<class T0, class Storage0, class Allocator0> 00187 void IlutPreconditioning<real, cplx, Allocator>:: 00188 FactorizeMatrix(const IVect& perm, 00189 Matrix<T0, General, Storage0, Allocator0>& mat, 00190 bool keep_matrix) 00191 { 00192 if (symmetric_algorithm) 00193 { 00194 cout << "Conversion to symmetric matrices not implemented." << endl; 00195 abort(); 00196 } 00197 else 00198 FactorizeUnsymMatrix(perm, mat, keep_matrix); 00199 } 00200 00201 00202 template<class real, class cplx, class Allocator> 00203 template<class T0, class Storage0, class Allocator0> 00204 void IlutPreconditioning<real, cplx, Allocator>:: 00205 FactorizeMatrix(const IVect& perm, 00206 Matrix<T0, Symmetric, Storage0, Allocator0>& mat, 00207 bool keep_matrix) 00208 { 00209 if (symmetric_algorithm) 00210 FactorizeSymMatrix(perm, mat, keep_matrix); 00211 else 00212 FactorizeUnsymMatrix(perm, mat, keep_matrix); 00213 } 00214 00215 00216 template<class real, class cplx, class Allocator> 00217 template<class MatrixSparse> 00218 void IlutPreconditioning<real, cplx, Allocator>:: 00219 FactorizeSymMatrix(const IVect& perm, MatrixSparse& mat, bool keep_matrix) 00220 { 00221 IVect inv_permutation; 00222 00223 // We convert matrix to symmetric format. 00224 Copy(mat, mat_sym); 00225 00226 // Old matrix is erased if needed. 00227 if (!keep_matrix) 00228 mat.Clear(); 00229 00230 // We keep permutation array in memory, and check it. 00231 int n = mat_sym.GetM(); 00232 if (perm.GetM() != n) 00233 { 00234 cout << "Numbering array should have the same size as matrix."; 00235 cout << endl; 00236 abort(); 00237 } 00238 00239 permutation_row.Reallocate(n); 00240 inv_permutation.Reallocate(n); 00241 inv_permutation.Fill(-1); 00242 for (int i = 0; i < n; i++) 00243 { 00244 permutation_row(i) = perm(i); 00245 inv_permutation(perm(i)) = i; 00246 } 00247 00248 for (int i = 0; i < n; i++) 00249 if (inv_permutation(i) == -1) 00250 { 00251 cout << "Error in the numbering array." << endl; 00252 abort(); 00253 } 00254 00255 // Matrix is permuted. 00256 ApplyInversePermutation(mat_sym, perm, perm); 00257 00258 // Temporary vector used for solving. 00259 xtmp.Reallocate(n); 00260 00261 // Factorization is performed. 00262 GetIlut(*this, mat_sym); 00263 } 00264 00265 00266 template<class real, class cplx, class Allocator> 00267 template<class MatrixSparse> 00268 void IlutPreconditioning<real, cplx, Allocator>:: 00269 FactorizeUnsymMatrix(const IVect& perm, MatrixSparse& mat, bool keep_matrix) 00270 { 00271 IVect inv_permutation; 00272 00273 // We convert matrix to unsymmetric format. 00274 Copy(mat, mat_unsym); 00275 00276 // Old matrix is erased if needed. 00277 if (!keep_matrix) 00278 mat.Clear(); 00279 00280 // We keep permutation array in memory, and check it. 00281 int n = mat_unsym.GetM(); 00282 if (perm.GetM() != n) 00283 { 00284 cout << "Numbering array should have the same size as matrix."; 00285 cout << endl; 00286 abort(); 00287 } 00288 00289 permutation_row.Reallocate(n); 00290 permutation_col.Reallocate(n); 00291 inv_permutation.Reallocate(n); 00292 inv_permutation.Fill(-1); 00293 for (int i = 0; i < n; i++) 00294 { 00295 permutation_row(i) = i; 00296 permutation_col(i) = i; 00297 inv_permutation(perm(i)) = i; 00298 } 00299 00300 for (int i = 0; i < n; i++) 00301 if (inv_permutation(i) == -1) 00302 { 00303 cout << "Error in the numbering array." << endl; 00304 abort(); 00305 } 00306 00307 IVect iperm = inv_permutation; 00308 00309 // Rows of matrix are permuted. 00310 ApplyInversePermutation(mat_unsym, perm, perm); 00311 00312 // Temporary vector used for solving. 00313 xtmp.Reallocate(n); 00314 00315 // Factorization is performed. 00316 // Columns are permuted during the factorization. 00317 inv_permutation.Fill(); 00318 GetIlut(*this, mat_unsym, permutation_col, inv_permutation); 00319 00320 // Combining permutations. 00321 IVect itmp = permutation_col; 00322 for (int i = 0; i < n; i++) 00323 permutation_col(i) = iperm(itmp(i)); 00324 00325 permutation_row = perm; 00326 } 00327 00328 00329 template<class real, class cplx, class Allocator> 00330 template<class Matrix1, class Vector1> 00331 void IlutPreconditioning<real, cplx, Allocator>:: 00332 Solve(const Matrix1& A, const Vector1& r, Vector1& z) 00333 { 00334 if (symmetric_algorithm) 00335 { 00336 for (int i = 0; i < r.GetM(); i++) 00337 xtmp(permutation_row(i)) = r(i); 00338 00339 SolveLU(mat_sym, xtmp); 00340 00341 for (int i = 0; i < r.GetM(); i++) 00342 z(i) = xtmp(permutation_row(i)); 00343 } 00344 else 00345 { 00346 for (int i = 0; i < r.GetM(); i++) 00347 xtmp(permutation_row(i)) = r(i); 00348 00349 SolveLU(mat_unsym, xtmp); 00350 00351 for (int i = 0; i < r.GetM(); i++) 00352 z(permutation_col(i)) = xtmp(i); 00353 } 00354 } 00355 00356 00357 template<class real, class cplx, class Allocator> 00358 template<class Matrix1, class Vector1> 00359 void IlutPreconditioning<real, cplx, Allocator>:: 00360 TransSolve(const Matrix1& A, const Vector1& r, Vector1& z) 00361 { 00362 if (symmetric_algorithm) 00363 Solve(A, r, z); 00364 else 00365 { 00366 for (int i = 0; i < r.GetM(); i++) 00367 xtmp(i) = r(permutation_col(i)); 00368 00369 SolveLU(SeldonTrans, mat_unsym, xtmp); 00370 00371 for (int i = 0; i < r.GetM(); i++) 00372 z(i) = xtmp(permutation_row(i)); 00373 } 00374 } 00375 00376 00377 template<class real, class cplx, class Allocator> 00378 template<class Vector1> 00379 void IlutPreconditioning<real, cplx, Allocator>::Solve(Vector1& z) 00380 { 00381 if (symmetric_algorithm) 00382 { 00383 for (int i = 0; i < z.GetM(); i++) 00384 xtmp(permutation_row(i)) = z(i); 00385 00386 SolveLU(mat_sym, xtmp); 00387 00388 for (int i = 0; i < z.GetM(); i++) 00389 z(i) = xtmp(permutation_row(i)); 00390 } 00391 else 00392 { 00393 for (int i = 0; i < z.GetM(); i++) 00394 xtmp(permutation_row(i)) = z(i); 00395 00396 SolveLU(mat_unsym, xtmp); 00397 00398 for (int i = 0; i < z.GetM(); i++) 00399 z(permutation_col(i)) = xtmp(i); 00400 } 00401 } 00402 00403 00404 template<class real, class cplx, class Allocator> 00405 template<class Vector1> 00406 void IlutPreconditioning<real, cplx, Allocator>::TransSolve(Vector1& z) 00407 { 00408 if (symmetric_algorithm) 00409 Solve(z); 00410 else 00411 { 00412 for (int i = 0; i < z.GetM(); i++) 00413 xtmp(i) = z(permutation_col(i)); 00414 00415 SolveLU(SeldonTrans, mat_unsym, xtmp); 00416 00417 for (int i = 0; i < z.GetM(); i++) 00418 z(i) = xtmp(permutation_row(i)); 00419 } 00420 } 00421 00422 00423 template<class real, class cplx, class Allocator> 00424 template<class TransStatus, class Vector1> 00425 void IlutPreconditioning<real, cplx, Allocator> 00426 ::Solve(const TransStatus& transA, Vector1& z) 00427 { 00428 if (transA.Trans()) 00429 TransSolve(z); 00430 else 00431 Solve(z); 00432 } 00433 00434 00435 template<class real, class cplx, class Storage, class Allocator> 00436 void qsplit_ilut(Vector<cplx, Storage, Allocator>& a, IVect& ind, int first, 00437 int n, int ncut, const real& abs_ncut) 00438 { 00439 //----------------------------------------------------------------------- 00440 // does a quick-sort split of a real array. 00441 // on input a(1:n). is a real array 00442 // on output a(1:n) is permuted such that its elements satisfy: 00443 // 00444 // abs(a(i)) .ge. abs(a(ncut)) for i .lt. ncut and 00445 // abs(a(i)) .le. abs(a(ncut)) for i .gt. ncut 00446 // 00447 // ind is an integer array which permuted in the same way as a 00448 //----------------------------------------------------------------------- 00449 int last = n-1; 00450 int ncut_ = ncut-1; 00451 int first_ = first; 00452 00453 if ((ncut_ < first_) || (ncut_ > last)) 00454 return; 00455 00456 cplx tmp; real abskey; 00457 int mid, itmp; 00458 bool test_loop = true; 00459 00460 // outer loop -- while mid .ne. ncut do 00461 while (test_loop) 00462 { 00463 mid = first_; 00464 abskey = abs(a(mid)); 00465 for (int j = (first_+1); j <= last; j++) 00466 { 00467 if (abs(a(j)) > abskey) 00468 { 00469 mid++; 00470 // Interchange. 00471 tmp = a(mid); 00472 itmp = ind(mid); 00473 a(mid) = a(j); 00474 ind(mid) = ind(j); 00475 a(j) = tmp; 00476 ind(j) = itmp; 00477 } 00478 } 00479 00480 // Interchange. 00481 tmp = a(mid); 00482 a(mid) = a(first_); 00483 a(first_) = tmp; 00484 00485 itmp = ind(mid); 00486 ind(mid) = ind(first_); 00487 ind(first_) = itmp; 00488 00489 // Test for while loop. 00490 if (mid == ncut_) 00491 return; 00492 00493 if (mid > ncut_) 00494 last = mid-1; 00495 else 00496 first_ = mid+1; 00497 } 00498 } 00499 00500 00502 template<class real, class cplx, class Allocator1, class Allocator2> 00503 void GetIlut(const IlutPreconditioning<real, cplx, Allocator1>& param, 00504 Matrix<cplx, General, ArrayRowSparse, Allocator2>& A, 00505 IVect& iperm, IVect& rperm) 00506 { 00507 int size_row; 00508 int n = A.GetN(); 00509 int type_factorization = param.GetFactorizationType(); 00510 int lfil = param.GetFillLevel(); 00511 real zero(0); 00512 real droptol = param.GetDroppingThreshold(); 00513 real alpha = param.GetDiagonalCoefficient(); 00514 bool variable_fill = false; 00515 bool standard_dropping = true; 00516 int additional_fill = param.GetAdditionalFillNumber(); 00517 int mbloc = param.GetPivotBlockInteger(); 00518 real permtol = param.GetPivotThreshold(); 00519 int print_level = param.GetPrintLevel(); 00520 00521 if (type_factorization == param.ILUT) 00522 standard_dropping = false; 00523 else if (type_factorization == param.ILU_D) 00524 standard_dropping = true; 00525 else if (type_factorization == param.ILUT_K) 00526 { 00527 variable_fill = true; // We use a variable lfil 00528 standard_dropping = false; 00529 } 00530 else if (type_factorization == param.ILU_0) 00531 { 00532 GetIlu0(A); 00533 return; 00534 } 00535 else if (type_factorization == param.MILU_0) 00536 { 00537 GetMilu0(A); 00538 return; 00539 } 00540 else if (type_factorization == param.ILU_K) 00541 { 00542 GetIluk(lfil, A); 00543 return; 00544 } 00545 00546 cplx fact, s, t; 00547 real tnorm; 00548 int length_lower, length_upper, jpos, jrow, i_row, j_col; 00549 int i, j, k, index_lu, length; 00550 00551 00552 if (lfil < 0) 00553 { 00554 cout << "Incorrect fill level." << endl; 00555 abort(); 00556 } 00557 00558 cplx czero, cone; 00559 SetComplexZero(czero); 00560 SetComplexOne(cone); 00561 typedef Vector<cplx, VectFull, Allocator2> VectCplx; 00562 VectCplx Row_Val(n); 00563 IVect Index(n), Row_Ind(n), Index_Diag(n); 00564 Row_Val.Fill(czero); 00565 Row_Ind.Fill(-1); 00566 Index_Diag.Fill(-1); 00567 00568 Index.Fill(-1); 00569 00570 bool element_dropped; cplx dropsum; 00571 00572 // main loop 00573 int new_percent = 0, old_percent = 0; 00574 for (i_row = 0; i_row < n; i_row++) 00575 { 00576 // Progress bar if print level is high enough. 00577 if (print_level > 0) 00578 { 00579 new_percent = int(double(i_row)/(n-1)*80); 00580 for (int percent = old_percent; percent < new_percent; percent++) 00581 { 00582 cout << "#"; cout.flush(); 00583 } 00584 00585 old_percent = new_percent; 00586 } 00587 00588 size_row = A.GetRowSize(i_row); 00589 tnorm = zero; 00590 00591 dropsum = czero; 00592 for (k = 0 ; k < size_row; k++) 00593 if (A.Value(i_row, k) != czero) 00594 tnorm += abs(A.Value(i_row, k)); 00595 00596 if (tnorm == zero) 00597 { 00598 cout << "Structurally singular matrix." << endl; 00599 cout << "Norm of row " << i_row << " is equal to 0." << endl; 00600 abort(); 00601 } 00602 00603 // tnorm is the sum of absolute value of coefficients of row i_row. 00604 tnorm /= real(size_row); 00605 if (variable_fill) 00606 lfil = size_row + additional_fill; 00607 00608 // Unpack L-part and U-part of row of A. 00609 length_upper = 1; 00610 length_lower = 0; 00611 Row_Ind(i_row) = i_row; 00612 Row_Val(i_row) = czero; 00613 Index(i_row) = i_row; 00614 00615 for (j = 0; j < size_row; j++) 00616 { 00617 k = rperm(A.Index(i_row, j)); 00618 00619 t = A.Value(i_row,j); 00620 if (k < i_row) 00621 { 00622 Row_Ind(length_lower) = k; 00623 Row_Val(length_lower) = t; 00624 Index(k) = length_lower; 00625 ++length_lower; 00626 } 00627 else if (k == i_row) 00628 { 00629 Row_Val(i_row) = t; 00630 } 00631 else 00632 { 00633 jpos = i_row + length_upper; 00634 Row_Ind(jpos) = k; 00635 Row_Val(jpos) = t; 00636 Index(k) = jpos; 00637 length_upper++; 00638 } 00639 } 00640 00641 j_col = 0; 00642 length = 0; 00643 00644 // Eliminates previous rows. 00645 while (j_col < length_lower) 00646 { 00647 // In order to do the elimination in the correct order, we must 00648 // select the smallest column index. 00649 jrow = Row_Ind(j_col); 00650 k = j_col; 00651 00652 // Determine smallest column index. 00653 for (j = j_col + 1; j < length_lower; j++) 00654 { 00655 if (Row_Ind(j) < jrow) 00656 { 00657 jrow = Row_Ind(j); 00658 k = j; 00659 } 00660 } 00661 00662 if (k != j_col) 00663 { 00664 // Exchanging column coefficients. 00665 j = Row_Ind(j_col); 00666 Row_Ind(j_col) = Row_Ind(k); 00667 Row_Ind(k) = j; 00668 00669 Index(jrow) = j_col; 00670 Index(j) = k; 00671 00672 s = Row_Val(j_col); 00673 Row_Val(j_col) = Row_Val(k); 00674 Row_Val(k) = s; 00675 } 00676 00677 // Zero out element in row. 00678 Index(jrow) = -1; 00679 00680 element_dropped = false; 00681 if (standard_dropping) 00682 if (abs(Row_Val(j_col)) <= droptol*tnorm) 00683 { 00684 dropsum += Row_Val(j_col); 00685 element_dropped = true; 00686 } 00687 00688 // Gets the multiplier for row to be eliminated (jrow). 00689 if (!element_dropped) 00690 { 00691 // first_index_upper points now on the diagonal coefficient. 00692 fact = Row_Val(j_col) * A.Value(jrow, Index_Diag(jrow)); 00693 00694 if (!standard_dropping) 00695 { 00696 if (abs(fact) <= droptol) 00697 element_dropped = true; 00698 } 00699 } 00700 00701 if (!element_dropped) 00702 { 00703 // Combines current row and row jrow. 00704 for (k = (Index_Diag(jrow)+1); k < A.GetRowSize(jrow); k++) 00705 { 00706 s = fact * A.Value(jrow,k); 00707 j = rperm(A.Index(jrow,k)); 00708 00709 jpos = Index(j); 00710 00711 if (j >= i_row) 00712 { 00713 00714 // Dealing with upper part. 00715 if (jpos == -1) 00716 { 00717 // This is a fill-in element. 00718 i = i_row + length_upper; 00719 Row_Ind(i) = j; 00720 Index(j) = i; 00721 Row_Val(i) = -s; 00722 ++length_upper; 00723 } 00724 else 00725 { 00726 // This is not a fill-in element. 00727 Row_Val(jpos) -= s; 00728 } 00729 } 00730 else 00731 { 00732 // Dealing with lower part. 00733 if (jpos == -1) 00734 { 00735 // this is a fill-in element 00736 Row_Ind(length_lower) = j; 00737 Index(j) = length_lower; 00738 Row_Val(length_lower) = -s; 00739 ++length_lower; 00740 } 00741 else 00742 { 00743 // This is not a fill-in element. 00744 Row_Val(jpos) -= s; 00745 } 00746 } 00747 } 00748 00749 // Stores this pivot element from left to right -- no danger 00750 // of overlap with the working elements in L (pivots). 00751 Row_Val(length) = fact; 00752 Row_Ind(length) = jrow; 00753 ++length; 00754 } 00755 00756 j_col++; 00757 } 00758 00759 // Resets double-pointer to zero (U-part). 00760 for (k = 0; k < length_upper; k++) 00761 Index(Row_Ind(i_row+k )) = -1; 00762 00763 // Updates L-matrix. 00764 if (!standard_dropping) 00765 { 00766 length_lower = length; 00767 length = min(length_lower,lfil); 00768 00769 // Sorts by quick-split. 00770 qsplit_ilut(Row_Val, Row_Ind, 0, length_lower, length,tnorm); 00771 } 00772 00773 size_row = length; 00774 A.ReallocateRow(i_row,size_row); 00775 00776 // store L-part 00777 index_lu = 0; 00778 for (k = 0 ; k < length ; k++) 00779 { 00780 A.Value(i_row,index_lu) = Row_Val(k); 00781 A.Index(i_row,index_lu) = iperm(Row_Ind(k)); 00782 ++index_lu; 00783 } 00784 00785 // Saves pointer to beginning of row i_row of U. 00786 Index_Diag(i_row) = index_lu; 00787 00788 // Updates. U-matrix -- first apply dropping strategy. 00789 length = 0; 00790 for (k = 1; k <= (length_upper-1); k++) 00791 { 00792 if (abs(Row_Val(i_row+k)) > droptol * tnorm) 00793 { 00794 ++length; 00795 Row_Val(i_row+length) = Row_Val(i_row+k); 00796 Row_Ind(i_row+length) = Row_Ind(i_row+k); 00797 } 00798 else 00799 dropsum += Row_Val(i_row+k); 00800 } 00801 00802 if (!standard_dropping) 00803 { 00804 length_upper = length + 1; 00805 length = min(length_upper,lfil); 00806 00807 qsplit_ilut(Row_Val, Row_Ind, i_row+1, 00808 i_row+length_upper, i_row+length+1, tnorm); 00809 } 00810 else 00811 length++; 00812 00813 // Determines next pivot. 00814 int imax = i_row; 00815 real xmax = abs(Row_Val(imax)); 00816 real xmax0 = xmax; 00817 int icut = i_row + mbloc - 1 - i_row%mbloc; 00818 for ( k = i_row + 1; k <= i_row + length - 1; k++) 00819 { 00820 tnorm = abs(Row_Val(k)); 00821 if ((tnorm > xmax) && (tnorm*permtol > xmax0) 00822 && (Row_Ind(k)<= icut)) 00823 { 00824 imax = k; 00825 xmax = tnorm; 00826 } 00827 } 00828 00829 // Exchanges Row_Val. 00830 s = Row_Val(i_row); 00831 Row_Val(i_row) = Row_Val(imax); 00832 Row_Val(imax) = s; 00833 00834 // Updates iperm and reverses iperm. 00835 j = Row_Ind(imax); 00836 i = iperm(i_row); 00837 iperm(i_row) = iperm(j); 00838 iperm(j) = i; 00839 00840 // Reverses iperm. 00841 rperm(iperm(i_row)) = i_row; 00842 rperm(iperm(j)) = j; 00843 00844 // Copies U-part in original coordinates. 00845 int index_diag = index_lu; 00846 A.ResizeRow(i_row, size_row+length); 00847 00848 for (k = i_row ; k <= i_row + length - 1; k++) 00849 { 00850 A.Index(i_row,index_lu) = iperm(Row_Ind(k)); 00851 A.Value(i_row,index_lu) = Row_Val(k); 00852 ++index_lu; 00853 } 00854 00855 // Stores inverse of diagonal element of u. 00856 if (standard_dropping) 00857 Row_Val(i_row) += alpha*dropsum; 00858 00859 if (Row_Val(i_row) == czero) 00860 Row_Val(i_row) = (droptol + 1e-4) * tnorm; 00861 00862 A.Value(i_row, index_diag) = cone / Row_Val(i_row); 00863 00864 } // end main loop 00865 00866 if (print_level > 0) 00867 cout << endl; 00868 00869 if (print_level > 0) 00870 cout << "The matrix takes " << 00871 int((A.GetDataSize()*(sizeof(cplx)+4))/(1024*1024)) << " MB" << endl; 00872 00873 for (i = 0; i < n; i++ ) 00874 for (j = 0; j < A.GetRowSize(i); j++) 00875 A.Index(i,j) = rperm(A.Index(i,j)); 00876 } 00877 00878 00879 template<class cplx, class Allocator> 00880 void GetIluk(int lfil, Matrix<cplx, General, ArrayRowSparse, Allocator>& A) 00881 { 00882 int n = A.GetM(); 00883 Vector<cplx, VectFull, CallocAlloc<cplx> > w; 00884 w.Reallocate(n+1); 00885 IVect jw(3*n), Index_Diag(n); 00886 Vector<IVect, VectFull, NewAlloc<IVect> > levs(n); 00887 00888 cplx czero, cone; 00889 SetComplexZero(czero); 00890 SetComplexOne(cone); 00891 00892 // Local variables 00893 cplx fact, s, t; 00894 int length_lower,length_upper, jpos, jrow, i_row, j_col; 00895 int i, j, k, index_lu, length; 00896 bool element_dropped; 00897 00898 int n2 = 2*n, jlev, k_, size_upper; 00899 jw.Fill(-1); 00900 00901 // Main loop. 00902 for (i_row = 0; i_row < n; i_row++) 00903 { 00904 int size_row = A.GetRowSize(i_row); 00905 00906 // Unpacks L-part and U-part of row of A in arrays w, jw. 00907 length_upper = 1; 00908 length_lower = 0; 00909 jw(i_row) = i_row; 00910 w(i_row) = 0.0; 00911 jw(n + i_row) = i_row; 00912 00913 for (j = 0; j < size_row; j++) 00914 { 00915 k = A.Index(i_row,j); 00916 t = A.Value(i_row,j); 00917 if (k < i_row) 00918 { 00919 jw(length_lower) = k; 00920 w(length_lower) = t; 00921 jw(n + k) = length_lower; 00922 jw(n2+length_lower) = -1; 00923 ++length_lower; 00924 } 00925 else if (k == i_row) 00926 { 00927 w(i_row) = t; 00928 jw(n2+length_lower) = -1; 00929 } 00930 else 00931 { 00932 jpos = i_row + length_upper; 00933 jw(jpos) = k; 00934 w(jpos) = t; 00935 jw(n + k) = jpos; 00936 length_upper++; 00937 } 00938 } 00939 00940 j_col = 0; 00941 length = 0; 00942 00943 00944 // Eliminates previous rows. 00945 while (j_col <length_lower) 00946 { 00947 // In order to do the elimination in the correct order, we must 00948 // select the smallest column index among jw(k); k = j_col + 1, 00949 // ..., length_lower. 00950 jrow = jw(j_col); 00951 k = j_col; 00952 00953 // Determines smallest column index. 00954 for (j = (j_col+1) ; j < length_lower; j++) 00955 { 00956 if (jw(j) < jrow) 00957 { 00958 jrow = jw(j); 00959 k = j; 00960 } 00961 } 00962 00963 if (k != j_col) 00964 { 00965 // Exchanges in jw. 00966 j = jw(j_col); 00967 jw(j_col) = jw(k); 00968 jw(k) = j; 00969 00970 // Exchanges in jw(n+ (pointers/ nonzero indicator). 00971 jw(n+jrow) = j_col; 00972 jw(n+j) = k; 00973 00974 // Exchanges in jw(n2+ (levels). 00975 j = jw(n2+j_col); 00976 jw(n2+j_col) = jw(n2+k); 00977 jw(n2+k) = j; 00978 00979 // Exchanges in w. 00980 s = w(j_col); 00981 w(j_col) = w(k); 00982 w(k) = s; 00983 } 00984 00985 // Zero out element in row by setting jw(n+jrow) to zero. 00986 jw(n + jrow) = -1; 00987 00988 element_dropped = false; 00989 00990 // Gets the multiplier for row to be eliminated (jrow). 00991 fact = w(j_col) * A.Value(jrow,Index_Diag(jrow)); 00992 00993 jlev = jw(n2+j_col) + 1; 00994 if (jlev > lfil) 00995 element_dropped = true; 00996 00997 if (!element_dropped) 00998 { 00999 // Combines current row and row jrow. 01000 k_ = 0; 01001 for (k = (Index_Diag(jrow)+1); k < A.GetRowSize(jrow) ; k++) 01002 { 01003 s = fact * A.Value(jrow,k); 01004 j = A.Index(jrow,k); 01005 01006 jpos = jw(n + j); 01007 if (j >= i_row) 01008 { 01009 // Dealing with upper part. 01010 if (jpos == -1) 01011 { 01012 // This is a fill-in element. 01013 i = i_row + length_upper; 01014 jw(i) = j; 01015 jw(n + j) = i; 01016 w(i) = -s; 01017 01018 jw(n2+i) = jlev + levs(jrow)(k_) + 1; 01019 ++length_upper; 01020 } 01021 else 01022 { 01023 // This is not a fill-in element. 01024 w(jpos) -= s; 01025 jw(n2+jpos) = min(jw(n2+jpos), 01026 jlev + levs(jrow)(k_)+1); 01027 } 01028 } 01029 else 01030 { 01031 // Dealing with lower part. 01032 if (jpos == -1) 01033 { 01034 // This is a fill-in element. 01035 jw(length_lower) = j; 01036 jw(n + j) = length_lower; 01037 w(length_lower) = -s; 01038 jw(n2+length_lower) = jlev + levs(jrow)(k_) + 1; 01039 ++length_lower; 01040 } 01041 else 01042 { 01043 // This is not a fill-in element. 01044 w(jpos) -= s; 01045 jw(n2+jpos) = min(jw(n2 + jpos), 01046 jlev + levs(jrow)(k_) + 1); 01047 } 01048 } 01049 01050 k_++; 01051 } 01052 01053 } 01054 01055 // Stores this pivot element from left to right -- no danger of 01056 // overlap with the working elements in L (pivots). 01057 w(j_col) = fact; 01058 jw(j_col) = jrow; 01059 01060 j_col++; 01061 } 01062 01063 // Resets double-pointer to zero (U-part). 01064 for (k = 0; k < length_upper; k++) 01065 jw(n + jw(i_row + k )) = -1; 01066 01067 // Updates L-matrix. 01068 size_row = 1; // we have the diagonal value. 01069 // Size of L-matrix. 01070 for (k = 0; k < length_lower; k++) 01071 if (jw(n2+k) < lfil) 01072 size_row++; 01073 01074 // Size of U-matrix. 01075 size_upper = 0; 01076 for (k = (i_row+1) ; k <= (i_row+length_upper-1) ; k++) 01077 if (jw(n2+k) < lfil) 01078 size_upper++; 01079 01080 size_row += size_upper; 01081 A.ReallocateRow(i_row,size_row); 01082 levs(i_row).Reallocate(size_upper); 01083 01084 index_lu = 0; 01085 for (k = 0; k < length_lower; k++) 01086 { 01087 if (jw(n2+k) < lfil) 01088 { 01089 A.Value(i_row,index_lu) = w(k); 01090 A.Index(i_row,index_lu) = jw(k); 01091 ++index_lu; 01092 } 01093 } 01094 01095 // Saves pointer to beginning of row i_row of U. 01096 Index_Diag(i_row) = index_lu; 01097 A.Value(i_row,index_lu) = cone / w(i_row); 01098 A.Index(i_row,index_lu++) = i_row; 01099 01100 for (k = (i_row+1) ; k <= (i_row+length_upper-1) ; k++) 01101 { 01102 if (jw(n2+k) < lfil) 01103 { 01104 A.Index(i_row,index_lu) = jw(k); 01105 A.Value(i_row,index_lu) = w(k); 01106 levs(i_row)(index_lu-Index_Diag(i_row)-1) = jw(n2+k); 01107 ++index_lu; 01108 } 01109 } 01110 } 01111 } 01112 01113 01114 template<class cplx, class Allocator> 01115 void GetIlu0(Matrix<cplx, General, ArrayRowSparse, Allocator>& A) 01116 { 01117 int j_col, jrow, jw, n = A.GetM(); 01118 IVect Index(n), ju(n); 01119 01120 cplx czero, cone; 01121 SetComplexZero(czero); 01122 SetComplexOne(cone); 01123 01124 // Initializes work vector to zero's. 01125 Index.Fill(-1); ju.Fill(-1); 01126 cplx tl; 01127 01128 // Main loop. 01129 for (int i_row = 0 ; i_row < n ; i_row++) 01130 { 01131 // Generating row number i_row of L and U. 01132 for (int j = 0 ; j < A.GetRowSize(i_row) ; j++ ) 01133 { 01134 j_col = A.Index(i_row, j); 01135 if (j_col == i_row) 01136 ju(i_row) = j; 01137 01138 Index(j_col) = j; 01139 } 01140 01141 int jm = ju(i_row)-1; 01142 01143 // Exit if diagonal element is reached. 01144 for (int j = 0; j <= jm; j++) 01145 { 01146 jrow = A.Index(i_row, j); 01147 tl = A.Value(i_row, j)*A.Value(jrow, ju(jrow)); 01148 A.Value(i_row, j) = tl; 01149 01150 // Performs linear combination. 01151 for ( j_col = (ju(jrow)+1); j_col < A.GetRowSize(jrow); j_col++) 01152 { 01153 jw = Index(A.Index(jrow,j_col)); 01154 if (jw != -1) 01155 A.Value(i_row, jw) -= tl*A.Value(jrow, j_col); 01156 } 01157 } 01158 01159 01160 // Inverts and stores diagonal element. 01161 if (A.Value(i_row, ju(i_row)) == czero) 01162 { 01163 cout << "Factorization fails because we found a null coefficient" 01164 << " on diagonal " << i_row << endl; 01165 abort(); 01166 } 01167 01168 A.Value(i_row,ju(i_row)) = cone / A.Value(i_row,ju(i_row)); 01169 01170 // Resets pointer Index to zero. 01171 Index(i_row) = -1; 01172 for (int i = 0; i < A.GetRowSize(i_row); i++) 01173 Index(A.Index(i_row, i)) = -1; 01174 } 01175 01176 } 01177 01178 01179 template<class cplx, class Allocator> 01180 void GetMilu0(Matrix<cplx, General, ArrayRowSparse, Allocator>& A) 01181 { 01182 int j_col, jrow, jw, n = A.GetM(); 01183 IVect Index(n), ju(n); 01184 01185 cplx czero, cone; 01186 SetComplexZero(czero); 01187 SetComplexOne(cone); 01188 01189 // Initializes work vector to zero's. 01190 Index.Fill(-1); ju.Fill(-1); 01191 cplx tl; 01192 01193 // Main loop. 01194 for (int i_row = 0 ; i_row < n ; i_row++) 01195 { 01196 // Generating row number i_row of L and U. 01197 for (int j = 0; j < A.GetRowSize(i_row); j++ ) 01198 { 01199 j_col = A.Index(i_row, j); 01200 if (j_col == i_row) 01201 ju(i_row) = j; 01202 01203 Index(j_col) = j; 01204 } 01205 01206 int jm = ju(i_row)-1; 01207 // Exit if diagonal element is reached. 01208 // s accumulates fill-in values. 01209 cplx s(0); 01210 for (int j = 0; j <= jm; j++) 01211 { 01212 jrow = A.Index(i_row, j); 01213 tl = A.Value(i_row, j)*A.Value(jrow, ju(jrow)); 01214 A.Value(i_row, j) = tl; 01215 01216 // Performs linear combination. 01217 for ( j_col = (ju(jrow)+1); j_col < A.GetRowSize(jrow); j_col++ ) 01218 { 01219 jw = Index(A.Index(jrow, j_col)); 01220 if (jw != -1) 01221 A.Value(i_row, jw) -= tl*A.Value(jrow, j_col); 01222 else 01223 s += tl*A.Value(jrow, j_col); 01224 } 01225 } 01226 01227 // Inverts and stores diagonal element. 01228 A.Value(i_row, ju(i_row)) -= s; 01229 if (A.Value(i_row, ju(i_row)) == czero) 01230 { 01231 cout << "Factorization fails because we found a null coefficient" 01232 << " on diagonal " << i_row << endl; 01233 abort(); 01234 } 01235 01236 A.Value(i_row, ju(i_row)) = cone /A.Value(i_row, ju(i_row)); 01237 01238 // Resets pointer Index to zero. 01239 Index(i_row) = -1; 01240 for (int i = 0; i < A.GetRowSize(i_row); i++) 01241 Index(A.Index(i_row, i)) = -1; 01242 } 01243 } 01244 01245 } 01246 01247 #define SELDON_FILE_ILUT_PRECONDITIONING_CXX 01248 #endif