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
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 #include "Ifpack_ConfigDefs.h"
00044 #include "Ifpack_Preconditioner.h"
00045 #include "Ifpack_ILUT.h"
00046 #include "Ifpack_Condest.h"
00047 #include "Ifpack_Utils.h"
00048 #include "Ifpack_HashTable.h"
00049 #include "Epetra_SerialComm.h"
00050 #include "Epetra_Comm.h"
00051 #include "Epetra_Map.h"
00052 #include "Epetra_RowMatrix.h"
00053 #include "Epetra_CrsMatrix.h"
00054 #include "Epetra_Vector.h"
00055 #include "Epetra_MultiVector.h"
00056 #include "Epetra_Util.h"
00057 #include "Teuchos_ParameterList.hpp"
00058 #include "Teuchos_RefCountPtr.hpp"
00059 #include <functional>
00060 #include <algorithm>
00061
00062 using namespace Teuchos;
00063
00064
00065
00066 Ifpack_ILUT::Ifpack_ILUT(const Epetra_RowMatrix* A) :
00067 A_(*A),
00068 Comm_(A->Comm()),
00069 Condest_(-1.0),
00070 Relax_(0.),
00071 Athresh_(0.0),
00072 Rthresh_(1.0),
00073 LevelOfFill_(1.0),
00074 DropTolerance_(1e-12),
00075 IsInitialized_(false),
00076 IsComputed_(false),
00077 UseTranspose_(false),
00078 NumMyRows_(-1),
00079 NumInitialize_(0),
00080 NumCompute_(0),
00081 NumApplyInverse_(0),
00082 InitializeTime_(0.0),
00083 ComputeTime_(0.0),
00084 ApplyInverseTime_(0.0),
00085 ComputeFlops_(0.0),
00086 ApplyInverseFlops_(0.0),
00087 Time_(Comm()),
00088 GlobalNonzeros_(0)
00089 {
00090
00091 }
00092
00093
00094 Ifpack_ILUT::~Ifpack_ILUT()
00095 {
00096 Destroy();
00097 }
00098
00099
00100 void Ifpack_ILUT::Destroy()
00101 {
00102 IsInitialized_ = false;
00103 IsComputed_ = false;
00104 }
00105
00106
00107 int Ifpack_ILUT::SetParameters(Teuchos::ParameterList& List)
00108 {
00109 try
00110 {
00111 LevelOfFill_ = List.get<double>("fact: ilut level-of-fill", LevelOfFill());
00112 if (LevelOfFill_ <= 0.0)
00113 IFPACK_CHK_ERR(-2);
00114
00115 Athresh_ = List.get<double>("fact: absolute threshold", Athresh_);
00116 Rthresh_ = List.get<double>("fact: relative threshold", Rthresh_);
00117 Relax_ = List.get<double>("fact: relax value", Relax_);
00118 DropTolerance_ = List.get<double>("fact: drop tolerance", DropTolerance_);
00119
00120 Label_ = "IFPACK ILUT (fill=" + Ifpack_toString(LevelOfFill())
00121 + ", relax=" + Ifpack_toString(RelaxValue())
00122 + ", athr=" + Ifpack_toString(AbsoluteThreshold())
00123 + ", rthr=" + Ifpack_toString(RelativeThreshold())
00124 + ", droptol=" + Ifpack_toString(DropTolerance())
00125 + ")";
00126 return(0);
00127 }
00128 catch (...)
00129 {
00130 cerr << "Caught an exception while parsing the parameter list" << endl;
00131 cerr << "This typically means that a parameter was set with the" << endl;
00132 cerr << "wrong type (for example, int instead of double). " << endl;
00133 cerr << "please check the documentation for the type required by each parameer." << endl;
00134 IFPACK_CHK_ERR(-1);
00135 }
00136
00137 return(0);
00138 }
00139
00140
00141 int Ifpack_ILUT::Initialize()
00142 {
00143
00144 Destroy();
00145
00146 Time_.ResetStartTime();
00147
00148
00149 if (Comm().NumProc() == 1 && Matrix().NumMyRows() != Matrix().NumMyCols())
00150 IFPACK_CHK_ERR(-2);
00151
00152 NumMyRows_ = Matrix().NumMyRows();
00153
00154
00155 IsInitialized_ = true;
00156 ++NumInitialize_;
00157 InitializeTime_ += Time_.ElapsedTime();
00158
00159 return(0);
00160 }
00161
00162
00163 class Ifpack_AbsComp
00164 {
00165 public:
00166 inline bool operator()(const double& x, const double& y)
00167 {
00168 return(IFPACK_ABS(x) > IFPACK_ABS(y));
00169 }
00170 };
00171
00172
00173
00174
00175
00176 template<typename int_type>
00177 int Ifpack_ILUT::TCompute()
00178 {
00179 int Length = A_.MaxNumEntries();
00180 std::vector<double> RowValuesL(Length);
00181 std::vector<int> RowIndicesU(Length);
00182 std::vector<int_type> RowIndicesU_LL;
00183 RowIndicesU_LL.resize(Length);
00184 std::vector<double> RowValuesU(Length);
00185
00186 int RowNnzU;
00187
00188 L_ = rcp(new Epetra_CrsMatrix(Copy, *SerialMap_, 0));
00189 U_ = rcp(new Epetra_CrsMatrix(Copy, *SerialMap_, 0));
00190
00191 if ((L_.get() == 0) || (U_.get() == 0))
00192 IFPACK_CHK_ERR(-5);
00193
00194
00195 IFPACK_CHK_ERR(A_.ExtractMyRowCopy(0,Length,RowNnzU,
00196 &RowValuesU[0],&RowIndicesU[0]));
00197
00198 bool distributed = (Comm().NumProc() > 1)?true:false;
00199
00200 if (distributed)
00201 {
00202 int count = 0;
00203 for (int i = 0 ;i < RowNnzU ; ++i)
00204 {
00205 if (RowIndicesU[i] < NumMyRows_){
00206 RowIndicesU[count] = RowIndicesU[i];
00207 RowValuesU[count] = RowValuesU[i];
00208 ++count;
00209 }
00210 else
00211 continue;
00212 }
00213 RowNnzU = count;
00214 }
00215
00216
00217 for (int i = 0 ;i < RowNnzU ; ++i) {
00218 if (RowIndicesU[i] == 0) {
00219 double& v = RowValuesU[i];
00220 v = AbsoluteThreshold() * EPETRA_SGN(v) + RelativeThreshold() * v;
00221 break;
00222 }
00223 }
00224
00225 std::copy(&(RowIndicesU[0]), &(RowIndicesU[0]) + RowNnzU, RowIndicesU_LL.begin());
00226 IFPACK_CHK_ERR(U_->InsertGlobalValues(0,RowNnzU,&(RowValuesU[0]),
00227 &(RowIndicesU_LL[0])));
00228
00229 RowValuesU[0] = 1.0;
00230 RowIndicesU[0] = 0;
00231 int_type RowIndicesU_0_templ = RowIndicesU[0];
00232 IFPACK_CHK_ERR(L_->InsertGlobalValues(0,1,&(RowValuesU[0]),
00233 &RowIndicesU_0_templ));
00234
00235 int max_keys = (int) 10 * A_.MaxNumEntries() * LevelOfFill() ;
00236 Ifpack_HashTable SingleRowU(max_keys , 1);
00237 Ifpack_HashTable SingleRowL(max_keys , 1);
00238
00239 int hash_size = SingleRowU.getRecommendedHashSize(max_keys) ;
00240 std::vector<int> keys; keys.reserve(hash_size * 10);
00241 std::vector<double> values; values.reserve(hash_size * 10);
00242 std::vector<double> AbsRow; AbsRow.reserve(hash_size * 10);
00243
00244
00245
00246
00247
00248 #ifdef IFPACK_FLOPCOUNTERS
00249 double this_proc_flops = 0.0;
00250 #endif
00251
00252 for (int row_i = 1 ; row_i < NumMyRows_ ; ++row_i)
00253 {
00254
00255 IFPACK_CHK_ERR(A_.ExtractMyRowCopy(row_i,Length,RowNnzU,
00256 &RowValuesU[0],&RowIndicesU[0]));
00257
00258 if (distributed)
00259 {
00260 int count = 0;
00261 for (int i = 0 ;i < RowNnzU ; ++i)
00262 {
00263 if (RowIndicesU[i] < NumMyRows_){
00264 RowIndicesU[count] = RowIndicesU[i];
00265 RowValuesU[count] = RowValuesU[i];
00266 ++count;
00267 }
00268 else
00269 continue;
00270 }
00271 RowNnzU = count;
00272 }
00273
00274 int NnzLower = 0;
00275 int NnzUpper = 0;
00276
00277 for (int i = 0 ;i < RowNnzU ; ++i) {
00278 if (RowIndicesU[i] < row_i)
00279 NnzLower++;
00280 else if (RowIndicesU[i] == row_i) {
00281
00282 NnzUpper++;
00283 double& v = RowValuesU[i];
00284 v = AbsoluteThreshold() * EPETRA_SGN(v) + RelativeThreshold() * v;
00285 }
00286 else
00287 NnzUpper++;
00288 }
00289
00290 int FillL = (int)(LevelOfFill() * NnzLower);
00291 int FillU = (int)(LevelOfFill() * NnzUpper);
00292 if (FillL == 0) FillL = 1;
00293 if (FillU == 0) FillU = 1;
00294
00295
00296 SingleRowU.reset();
00297
00298 for (int i = 0 ; i < RowNnzU ; ++i) {
00299 SingleRowU.set(RowIndicesU[i], RowValuesU[i]);
00300 }
00301
00302
00303 SingleRowL.reset();
00304
00305 int start_col = NumMyRows_;
00306 for (int i = 0 ; i < RowNnzU ; ++i)
00307 start_col = EPETRA_MIN(start_col, RowIndicesU[i]);
00308
00309 #ifdef IFPACK_FLOPCOUNTERS
00310 int flops = 0;
00311 #endif
00312
00313 for (int col_k = start_col ; col_k < row_i ; ++col_k) {
00314
00315 int_type* ColIndicesK;
00316 double* ColValuesK;
00317 int ColNnzK;
00318
00319 double xxx = SingleRowU.get(col_k);
00320
00321
00322 if (IFPACK_ABS(xxx) > DropTolerance()) {
00323 IFPACK_CHK_ERR(U_->ExtractGlobalRowView(col_k, ColNnzK, ColValuesK,
00324 ColIndicesK));
00325
00326
00327 double DiagonalValueK = 0.0;
00328 for (int i = 0 ; i < ColNnzK ; ++i) {
00329 if (ColIndicesK[i] == col_k) {
00330 DiagonalValueK = ColValuesK[i];
00331 break;
00332 }
00333 }
00334
00335
00336 if (DiagonalValueK == 0.0)
00337 DiagonalValueK = AbsoluteThreshold();
00338
00339 double yyy = xxx / DiagonalValueK ;
00340 SingleRowL.set(col_k, yyy);
00341 #ifdef IFPACK_FLOPCOUNTERS
00342 ++flops;
00343 #endif
00344
00345 for (int j = 0 ; yyy != 0.0 && j < ColNnzK ; ++j)
00346 {
00347 int_type col_j = ColIndicesK[j];
00348
00349 if (col_j < col_k) continue;
00350
00351 SingleRowU.set(col_j, -yyy * ColValuesK[j], true);
00352 #ifdef IFPACK_FLOPCOUNTERS
00353 flops += 2;
00354 #endif
00355 }
00356 }
00357 }
00358
00359 #ifdef IFPACK_FLOPCOUNTERS
00360 this_proc_flops += 1.0 * flops;
00361 #endif
00362
00363 double cutoff = DropTolerance();
00364 double DiscardedElements = 0.0;
00365 int count;
00366
00367
00368 count = 0;
00369 int sizeL = SingleRowL.getNumEntries();
00370 keys.resize(sizeL);
00371 values.resize(sizeL);
00372
00373 AbsRow.resize(sizeL);
00374
00375 SingleRowL.arrayify(
00376 keys.size() ? &keys[0] : 0,
00377 values.size() ? &values[0] : 0
00378 );
00379 for (int i = 0; i < sizeL; ++i)
00380 if (IFPACK_ABS(values[i]) > DropTolerance()) {
00381 AbsRow[count++] = IFPACK_ABS(values[i]);
00382 }
00383
00384 if (count > FillL) {
00385 nth_element(AbsRow.begin(), AbsRow.begin() + FillL, AbsRow.begin() + count,
00386 std::greater<double>());
00387 cutoff = AbsRow[FillL];
00388 }
00389
00390 for (int i = 0; i < sizeL; ++i) {
00391 if (IFPACK_ABS(values[i]) >= cutoff) {
00392 int_type key_templ = keys[i];
00393 IFPACK_CHK_ERR(L_->InsertGlobalValues(row_i,1, &values[i], &key_templ));
00394 }
00395 else
00396 DiscardedElements += values[i];
00397 }
00398
00399
00400
00401 double dtmp = 1.0;
00402 int_type row_i_templ = row_i;
00403 IFPACK_CHK_ERR(L_->InsertGlobalValues(row_i,1, &dtmp, &row_i_templ));
00404
00405
00406 count = 0;
00407 int sizeU = SingleRowU.getNumEntries();
00408 AbsRow.resize(sizeU + 1);
00409 keys.resize(sizeU + 1);
00410 values.resize(sizeU + 1);
00411
00412 SingleRowU.arrayify(&keys[0], &values[0]);
00413
00414 for (int i = 0; i < sizeU; ++i)
00415 if (keys[i] >= row_i && IFPACK_ABS(values[i]) > DropTolerance())
00416 {
00417 AbsRow[count++] = IFPACK_ABS(values[i]);
00418 }
00419
00420 if (count > FillU) {
00421 nth_element(AbsRow.begin(), AbsRow.begin() + FillU, AbsRow.begin() + count,
00422 std::greater<double>());
00423 cutoff = AbsRow[FillU];
00424 }
00425
00426
00427 for (int i = 0; i < sizeU; ++i)
00428 {
00429 int col = keys[i];
00430 double val = values[i];
00431
00432 if (col >= row_i) {
00433 if (IFPACK_ABS(val) >= cutoff || row_i == col) {
00434 int_type col_templ = col;
00435 IFPACK_CHK_ERR(U_->InsertGlobalValues(row_i,1, &val, &col_templ));
00436 }
00437 else
00438 DiscardedElements += val;
00439 }
00440 }
00441
00442
00443 if (RelaxValue() != 0.0) {
00444 DiscardedElements *= RelaxValue();
00445 IFPACK_CHK_ERR(U_->InsertGlobalValues(row_i,1, &DiscardedElements,
00446 &row_i_templ));
00447 }
00448 }
00449
00450 #ifdef IFPACK_FLOPCOUNTERS
00451 double tf;
00452 Comm().SumAll(&this_proc_flops, &tf, 1);
00453 ComputeFlops_ += tf;
00454 #endif
00455
00456 IFPACK_CHK_ERR(L_->FillComplete());
00457 IFPACK_CHK_ERR(U_->FillComplete());
00458
00459 #if 0
00460
00461 Epetra_Vector LHS(A_.RowMatrixRowMap());
00462 Epetra_Vector RHS1(A_.RowMatrixRowMap());
00463 Epetra_Vector RHS2(A_.RowMatrixRowMap());
00464 Epetra_Vector RHS3(A_.RowMatrixRowMap());
00465 LHS.Random();
00466
00467 cout << "A = " << A_.NumGlobalNonzeros() << endl;
00468 cout << "L = " << L_->NumGlobalNonzeros() << endl;
00469 cout << "U = " << U_->NumGlobalNonzeros() << endl;
00470
00471 A_.Multiply(false,LHS,RHS1);
00472 U_->Multiply(false,LHS,RHS2);
00473 L_->Multiply(false,RHS2,RHS3);
00474
00475 RHS1.Update(-1.0, RHS3, 1.0);
00476 double Norm;
00477 RHS1.Norm2(&Norm);
00478 #endif
00479
00480 long long MyNonzeros = L_->NumGlobalNonzeros64() + U_->NumGlobalNonzeros64();
00481 Comm().SumAll(&MyNonzeros, &GlobalNonzeros_, 1);
00482
00483 IsComputed_ = true;
00484
00485 ++NumCompute_;
00486 ComputeTime_ += Time_.ElapsedTime();
00487
00488 return(0);
00489
00490 }
00491
00492 int Ifpack_ILUT::Compute() {
00493 if (!IsInitialized())
00494 IFPACK_CHK_ERR(Initialize());
00495
00496 Time_.ResetStartTime();
00497 IsComputed_ = false;
00498
00499 NumMyRows_ = A_.NumMyRows();
00500 bool distributed = (Comm().NumProc() > 1)?true:false;
00501
00502 #if !defined(EPETRA_NO_32BIT_GLOBAL_INDICES) || !defined(EPETRA_NO_64BIT_GLOBAL_INDICES)
00503 if (distributed)
00504 {
00505 SerialComm_ = rcp(new Epetra_SerialComm);
00506 SerialMap_ = rcp(new Epetra_Map(NumMyRows_, 0, *SerialComm_));
00507 assert (SerialComm_.get() != 0);
00508 assert (SerialMap_.get() != 0);
00509 }
00510 else
00511 SerialMap_ = rcp(const_cast<Epetra_Map*>(&A_.RowMatrixRowMap()), false);
00512 #endif
00513
00514 int ret_val = 1;
00515
00516 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00517 if(SerialMap_->GlobalIndicesInt()) {
00518 ret_val = TCompute<int>();
00519 }
00520 else
00521 #endif
00522 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00523 if(SerialMap_->GlobalIndicesLongLong()) {
00524 ret_val = TCompute<long long>();
00525 }
00526 else
00527 #endif
00528 throw "Ifpack_ILUT::Compute: GlobalIndices type unknown";
00529
00530 return ret_val;
00531 }
00532
00533
00534 int Ifpack_ILUT::ApplyInverse(const Epetra_MultiVector& X,
00535 Epetra_MultiVector& Y) const
00536 {
00537 if (!IsComputed())
00538 IFPACK_CHK_ERR(-2);
00539
00540 if (X.NumVectors() != Y.NumVectors())
00541 IFPACK_CHK_ERR(-3);
00542
00543 Time_.ResetStartTime();
00544
00545
00546
00547
00548
00549
00550
00551 Teuchos::RefCountPtr<const Epetra_MultiVector> Xcopy;
00552 if (X.Pointers()[0] == Y.Pointers()[0])
00553 Xcopy = Teuchos::rcp( new Epetra_MultiVector(X) );
00554 else
00555 Xcopy = Teuchos::rcp( &X, false );
00556
00557 if (!UseTranspose_)
00558 {
00559
00560 IFPACK_CHK_ERR(L_->Solve(false,false,false,*Xcopy,Y));
00561 IFPACK_CHK_ERR(U_->Solve(true,false,false,Y,Y));
00562 }
00563 else
00564 {
00565
00566 IFPACK_CHK_ERR(U_->Solve(true,true,false,*Xcopy,Y));
00567 IFPACK_CHK_ERR(L_->Solve(false,true,false,Y,Y));
00568 }
00569
00570 ++NumApplyInverse_;
00571 #ifdef IFPACK_FLOPCOUNTERS
00572 ApplyInverseFlops_ += X.NumVectors() * 2 * GlobalNonzeros_;
00573 #endif
00574 ApplyInverseTime_ += Time_.ElapsedTime();
00575
00576 return(0);
00577
00578 }
00579
00580
00581 int Ifpack_ILUT::Apply(const Epetra_MultiVector& X,
00582 Epetra_MultiVector& Y) const
00583 {
00584 return(-98);
00585 }
00586
00587
00588 double Ifpack_ILUT::Condest(const Ifpack_CondestType CT,
00589 const int MaxIters, const double Tol,
00590 Epetra_RowMatrix* Matrix_in)
00591 {
00592 if (!IsComputed())
00593 return(-1.0);
00594
00595
00596 if (Condest_ == -1.0)
00597 Condest_ = Ifpack_Condest(*this, CT, MaxIters, Tol, Matrix_in);
00598
00599 return(Condest_);
00600 }
00601
00602
00603 std::ostream&
00604 Ifpack_ILUT::Print(std::ostream& os) const
00605 {
00606 if (!Comm().MyPID()) {
00607 os << endl;
00608 os << "================================================================================" << endl;
00609 os << "Ifpack_ILUT: " << Label() << endl << endl;
00610 os << "Level-of-fill = " << LevelOfFill() << endl;
00611 os << "Absolute threshold = " << AbsoluteThreshold() << endl;
00612 os << "Relative threshold = " << RelativeThreshold() << endl;
00613 os << "Relax value = " << RelaxValue() << endl;
00614 os << "Condition number estimate = " << Condest() << endl;
00615 os << "Global number of rows = " << A_.NumGlobalRows64() << endl;
00616 if (IsComputed_) {
00617 os << "Number of nonzeros in A = " << A_.NumGlobalNonzeros64() << endl;
00618 os << "Number of nonzeros in L + U = " << NumGlobalNonzeros64()
00619 << " ( = " << 100.0 * NumGlobalNonzeros64() / A_.NumGlobalNonzeros64()
00620 << " % of A)" << endl;
00621 os << "nonzeros / rows = "
00622 << 1.0 * NumGlobalNonzeros64() / U_->NumGlobalRows64() << endl;
00623 }
00624 os << endl;
00625 os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
00626 os << "----- ------- -------------- ------------ --------" << endl;
00627 os << "Initialize() " << std::setw(5) << NumInitialize()
00628 << " " << std::setw(15) << InitializeTime()
00629 << " 0.0 0.0" << endl;
00630 os << "Compute() " << std::setw(5) << NumCompute()
00631 << " " << std::setw(15) << ComputeTime()
00632 << " " << std::setw(15) << 1.0e-6 * ComputeFlops();
00633 if (ComputeTime() != 0.0)
00634 os << " " << std::setw(15) << 1.0e-6 * ComputeFlops() / ComputeTime() << endl;
00635 else
00636 os << " " << std::setw(15) << 0.0 << endl;
00637 os << "ApplyInverse() " << std::setw(5) << NumApplyInverse()
00638 << " " << std::setw(15) << ApplyInverseTime()
00639 << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops();
00640 if (ApplyInverseTime() != 0.0)
00641 os << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops() / ApplyInverseTime() << endl;
00642 else
00643 os << " " << std::setw(15) << 0.0 << endl;
00644 os << "================================================================================" << endl;
00645 os << endl;
00646 }
00647
00648 return(os);
00649 }