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_CondestType.h"
00045 #include "Ifpack_ILU.h"
00046 #include "Epetra_ConfigDefs.h"
00047 #include "Epetra_Comm.h"
00048 #include "Epetra_Map.h"
00049 #include "Epetra_RowMatrix.h"
00050 #include "Epetra_Vector.h"
00051 #include "Epetra_MultiVector.h"
00052 #include "Epetra_CrsGraph.h"
00053 #include "Epetra_CrsMatrix.h"
00054 #include "Teuchos_ParameterList.hpp"
00055 #include "Teuchos_RefCountPtr.hpp"
00056
00057 using Teuchos::RefCountPtr;
00058 using Teuchos::rcp;
00059
00060 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00061 # include "Teuchos_TimeMonitor.hpp"
00062 #endif
00063
00064
00065 Ifpack_ILU::Ifpack_ILU(Epetra_RowMatrix* Matrix_in) :
00066 A_(rcp(Matrix_in,false)),
00067 Comm_(Matrix_in->Comm()),
00068 UseTranspose_(false),
00069 NumMyDiagonals_(0),
00070 RelaxValue_(0.0),
00071 Athresh_(0.0),
00072 Rthresh_(1.0),
00073 Condest_(-1.0),
00074 LevelOfFill_(0),
00075 IsInitialized_(false),
00076 IsComputed_(false),
00077 NumInitialize_(0),
00078 NumCompute_(0),
00079 NumApplyInverse_(0),
00080 InitializeTime_(0.0),
00081 ComputeTime_(0.0),
00082 ApplyInverseTime_(0.0),
00083 ComputeFlops_(0.0),
00084 ApplyInverseFlops_(0.0),
00085 Time_(Comm())
00086 {
00087 Teuchos::ParameterList List;
00088 SetParameters(List);
00089 }
00090
00091
00092 void Ifpack_ILU::Destroy()
00093 {
00094
00095 U_DomainMap_ = 0;
00096 L_RangeMap_ = 0;
00097 }
00098
00099
00100 int Ifpack_ILU::SetParameters(Teuchos::ParameterList& List)
00101 {
00102 RelaxValue_ = List.get("fact: relax value", RelaxValue_);
00103 Athresh_ = List.get("fact: absolute threshold", Athresh_);
00104 Rthresh_ = List.get("fact: relative threshold", Rthresh_);
00105 LevelOfFill_ = List.get("fact: level-of-fill", LevelOfFill_);
00106
00107
00108 sprintf(Label_, "IFPACK ILU (fill=%d, relax=%f, athr=%f, rthr=%f)",
00109 LevelOfFill(), RelaxValue(), AbsoluteThreshold(),
00110 RelativeThreshold());
00111 return(0);
00112 }
00113
00114
00115 int Ifpack_ILU::ComputeSetup()
00116 {
00117
00118 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00119 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::ComputeSetup");
00120 #endif
00121
00122 L_ = rcp(new Epetra_CrsMatrix(Copy, Graph().L_Graph()));
00123 U_ = rcp(new Epetra_CrsMatrix(Copy, Graph().U_Graph()));
00124 D_ = rcp(new Epetra_Vector(Graph().L_Graph().RowMap()));
00125 if ((L_.get() == 0) || (U_.get() == 0) || (D_.get() == 0))
00126 IFPACK_CHK_ERR(-5);
00127
00128
00129 int MaxNumEntries = Matrix().MaxNumEntries();
00130
00131
00132 U_DomainMap_ = &(Matrix().OperatorDomainMap());
00133 L_RangeMap_ = &(Matrix().OperatorRangeMap());
00134
00135
00136 int ierr = 0;
00137 int i, j;
00138 int NumIn, NumL, NumU;
00139 bool DiagFound;
00140 int NumNonzeroDiags = 0;
00141
00142 std::vector<int> InI(MaxNumEntries);
00143 std::vector<int> LI(MaxNumEntries);
00144 std::vector<int> UI(MaxNumEntries);
00145 std::vector<double> InV(MaxNumEntries);
00146 std::vector<double> LV(MaxNumEntries);
00147 std::vector<double> UV(MaxNumEntries);
00148
00149 bool ReplaceValues = (L_->StaticGraph() || L_->IndicesAreLocal());
00150
00151 if (ReplaceValues) {
00152 L_->PutScalar(0.0);
00153 U_->PutScalar(0.0);
00154 }
00155
00156 D_->PutScalar(0.0);
00157 double *DV;
00158 IFPACK_CHK_ERR(D_->ExtractView(&DV));
00159
00160
00161
00162 for (i=0; i< NumMyRows(); i++) {
00163
00164 IFPACK_CHK_ERR(Matrix().ExtractMyRowCopy(i, MaxNumEntries, NumIn, &InV[0], &InI[0]));
00165
00166
00167
00168 NumL = 0;
00169 NumU = 0;
00170 DiagFound = false;
00171
00172 for (j=0; j< NumIn; j++) {
00173 int k = InI[j];
00174
00175 if (k==i) {
00176 DiagFound = true;
00177 DV[i] += Rthresh_ * InV[j] + EPETRA_SGN(InV[j]) * Athresh_;
00178 }
00179
00180 else if (k < 0) {IFPACK_CHK_ERR(-4);}
00181
00182 else if (k < i) {
00183 LI[NumL] = k;
00184 LV[NumL] = InV[j];
00185 NumL++;
00186 }
00187 else if (k<NumMyRows()) {
00188 UI[NumU] = k;
00189 UV[NumU] = InV[j];
00190 NumU++;
00191 }
00192 }
00193
00194
00195
00196 if (DiagFound) NumNonzeroDiags++;
00197 else DV[i] = Athresh_;
00198
00199 if (NumL) {
00200 if (ReplaceValues) {
00201 (L_->ReplaceMyValues(i, NumL, &LV[0], &LI[0]));
00202 }
00203 else {
00204 IFPACK_CHK_ERR(L_->InsertMyValues(i, NumL, &LV[0], &LI[0]));
00205 }
00206 }
00207
00208 if (NumU) {
00209 if (ReplaceValues) {
00210 (U_->ReplaceMyValues(i, NumU, &UV[0], &UI[0]));
00211 }
00212 else {
00213 IFPACK_CHK_ERR(U_->InsertMyValues(i, NumU, &UV[0], &UI[0]));
00214 }
00215 }
00216
00217 }
00218
00219 if (!ReplaceValues) {
00220
00221
00222
00223
00224 IFPACK_CHK_ERR(L_->FillComplete((L_->RowMatrixColMap()), *L_RangeMap_));
00225 IFPACK_CHK_ERR(U_->FillComplete(*U_DomainMap_, U_->RowMatrixRowMap()));
00226 }
00227
00228
00229 int TotalNonzeroDiags = 0;
00230 IFPACK_CHK_ERR(Graph().L_Graph().RowMap().Comm().SumAll(&NumNonzeroDiags, &TotalNonzeroDiags, 1));
00231 NumMyDiagonals_ = NumNonzeroDiags;
00232 if (NumNonzeroDiags != NumMyRows()) ierr = 1;
00233
00234 IFPACK_CHK_ERR(ierr);
00235 return(ierr);
00236 }
00237
00238
00239 int Ifpack_ILU::Initialize()
00240 {
00241
00242 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00243 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::Initialize");
00244 #endif
00245
00246 Time_.ResetStartTime();
00247 IsInitialized_ = false;
00248
00249
00250 Destroy();
00251
00252 Epetra_CrsMatrix* CrsMatrix;
00253 CrsMatrix = dynamic_cast<Epetra_CrsMatrix*>(&*A_);
00254 if (CrsMatrix == 0) {
00255
00256
00257
00258
00259 int size = A_->MaxNumEntries();
00260 CrsGraph_ = rcp(new Epetra_CrsGraph(Copy,A_->RowMatrixRowMap(), size));
00261 if (CrsGraph_.get() == 0)
00262 IFPACK_CHK_ERR(-5);
00263
00264 std::vector<double> Values(size);
00265
00266 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00267 if(A_->RowMatrixRowMap().GlobalIndicesInt()) {
00268 std::vector<int> Indices(size);
00269
00270
00271 for (int i = 0 ; i < A_->NumMyRows() ; ++i) {
00272 int NumEntries;
00273 int GlobalRow = A_->RowMatrixRowMap().GID(i);
00274 IFPACK_CHK_ERR(A_->ExtractMyRowCopy(i, size, NumEntries,
00275 &Values[0], &Indices[0]));
00276
00277 for (int j = 0 ; j < NumEntries ; ++j) {
00278 Indices[j] = A_->RowMatrixColMap().GID(Indices[j]);
00279 }
00280 IFPACK_CHK_ERR(CrsGraph_->InsertGlobalIndices(GlobalRow,NumEntries,
00281 &Indices[0]));
00282 }
00283 }
00284 else
00285 #endif
00286 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00287 if(A_->RowMatrixRowMap().GlobalIndicesLongLong()) {
00288 std::vector<int> Indices_local(size);
00289 std::vector<long long> Indices(size);
00290
00291
00292 for (int i = 0 ; i < A_->NumMyRows() ; ++i) {
00293 int NumEntries;
00294 long long GlobalRow = A_->RowMatrixRowMap().GID64(i);
00295 IFPACK_CHK_ERR(A_->ExtractMyRowCopy(i, size, NumEntries,
00296 &Values[0], &Indices_local[0]));
00297
00298 for (int j = 0 ; j < NumEntries ; ++j) {
00299 Indices[j] = A_->RowMatrixColMap().GID64(Indices_local[j]);
00300 }
00301 IFPACK_CHK_ERR(CrsGraph_->InsertGlobalIndices(GlobalRow,NumEntries,
00302 &Indices[0]));
00303 }
00304 }
00305 else
00306 #endif
00307 throw "Ifpack_ILU::Initialize: GlobalIndices type unknown";
00308
00309 IFPACK_CHK_ERR(CrsGraph_->FillComplete(A_->RowMatrixRowMap(),
00310 A_->RowMatrixRowMap()));
00311
00312
00313
00314 Graph_ = rcp(new Ifpack_IlukGraph(*CrsGraph_, LevelOfFill_, 0));
00315
00316 }
00317 else {
00318
00319 Graph_ = rcp(new Ifpack_IlukGraph(CrsMatrix->Graph(), LevelOfFill_, 0));
00320 }
00321
00322 if (Graph_.get() == 0)
00323 IFPACK_CHK_ERR(-5);
00324 IFPACK_CHK_ERR(Graph_->ConstructFilledGraph());
00325
00326 IsInitialized_ = true;
00327 NumInitialize_++;
00328 InitializeTime_ += Time_.ElapsedTime();
00329
00330 return(0);
00331 }
00332
00333
00334 int Ifpack_ILU::Compute()
00335 {
00336
00337 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00338 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::Compute");
00339 #endif
00340
00341 if (!IsInitialized())
00342 IFPACK_CHK_ERR(Initialize());
00343
00344 Time_.ResetStartTime();
00345 IsComputed_ = false;
00346
00347
00348 IFPACK_CHK_ERR(ComputeSetup());
00349
00350
00351
00352
00353 double MinDiagonalValue = Epetra_MinDouble;
00354 double MaxDiagonalValue = 1.0/MinDiagonalValue;
00355
00356 int ierr = 0;
00357 int i, j, k;
00358 int *LI, *UI;
00359 double *LV, *UV;
00360 int NumIn, NumL, NumU;
00361
00362
00363 int MaxNumEntries = L_->MaxNumEntries() + U_->MaxNumEntries() + 1;
00364
00365 std::vector<int> InI(MaxNumEntries+1);
00366 std::vector<double> InV(MaxNumEntries+1);
00367 std::vector<int> colflag(NumMyCols());
00368
00369 double *DV;
00370 ierr = D_->ExtractView(&DV);
00371
00372 #ifdef IFPACK_FLOPCOUNTERS
00373 int current_madds = 0;
00374 #endif
00375
00376
00377
00378
00379
00380
00381 int NumUU;
00382 int * UUI;
00383 double * UUV;
00384 for (j = 0; j < NumMyCols(); ++j) colflag[j] = - 1;
00385
00386 for (i = 0; i < NumMyRows(); ++i) {
00387
00388
00389
00390 NumIn = MaxNumEntries;
00391 IFPACK_CHK_ERR(L_->ExtractMyRowCopy(i, NumIn, NumL, &InV[0], &InI[0]));
00392 LV = &InV[0];
00393 LI = &InI[0];
00394
00395 InV[NumL] = DV[i];
00396 InI[NumL] = i;
00397
00398 IFPACK_CHK_ERR(U_->ExtractMyRowCopy(i, NumIn-NumL-1, NumU, &InV[NumL+1], &InI[NumL+1]));
00399 NumIn = NumL+NumU+1;
00400 UV = &InV[NumL+1];
00401 UI = &InI[NumL+1];
00402
00403
00404 for (j=0; j<NumIn; j++) colflag[InI[j]] = j;
00405
00406 double diagmod = 0.0;
00407
00408 for (int jj=0; jj<NumL; jj++) {
00409 j = InI[jj];
00410 double multiplier = InV[jj];
00411
00412 InV[jj] *= DV[j];
00413
00414 IFPACK_CHK_ERR(U_->ExtractMyRowView(j, NumUU, UUV, UUI));
00415
00416 if (RelaxValue_==0.0) {
00417 for (k=0; k<NumUU; k++) {
00418 int kk = colflag[UUI[k]];
00419 if (kk>-1) {
00420 InV[kk] -= multiplier*UUV[k];
00421 #ifdef IFPACK_FLOPCOUNTERS
00422 current_madds++;
00423 #endif
00424 }
00425 }
00426 }
00427 else {
00428 for (k=0; k<NumUU; k++) {
00429 int kk = colflag[UUI[k]];
00430 if (kk>-1) InV[kk] -= multiplier*UUV[k];
00431 else diagmod -= multiplier*UUV[k];
00432 #ifdef IFPACK_FLOPCOUNTERS
00433 current_madds++;
00434 #endif
00435 }
00436 }
00437 }
00438 if (NumL) {
00439 IFPACK_CHK_ERR(L_->ReplaceMyValues(i, NumL, LV, LI));
00440 }
00441
00442 DV[i] = InV[NumL];
00443
00444 if (RelaxValue_!=0.0) {
00445 DV[i] += RelaxValue_*diagmod;
00446
00447 }
00448
00449 if (fabs(DV[i]) > MaxDiagonalValue) {
00450 if (DV[i] < 0) DV[i] = - MinDiagonalValue;
00451 else DV[i] = MinDiagonalValue;
00452 }
00453 else
00454 DV[i] = 1.0/DV[i];
00455
00456 for (j=0; j<NumU; j++) UV[j] *= DV[i];
00457
00458 if (NumU) {
00459 IFPACK_CHK_ERR(U_->ReplaceMyValues(i, NumU, UV, UI));
00460 }
00461
00462
00463 for (j=0; j<NumIn; j++) colflag[InI[j]] = -1;
00464 }
00465
00466
00467
00468 if (!L_->LowerTriangular())
00469 IFPACK_CHK_ERR(-4);
00470 if (!U_->UpperTriangular())
00471 IFPACK_CHK_ERR(-4);
00472
00473 #ifdef IFPACK_FLOPCOUNTERS
00474
00475
00476 double current_flops = 2 * current_madds;
00477 double total_flops = 0;
00478
00479 IFPACK_CHK_ERR(Graph().L_Graph().RowMap().Comm().SumAll(¤t_flops, &total_flops, 1));
00480
00481
00482 total_flops += (double) L_->NumGlobalNonzeros();
00483 total_flops += (double) D_->GlobalLength();
00484 if (RelaxValue_!=0.0) total_flops += 2 * (double)D_->GlobalLength();
00485
00486
00487 ComputeFlops_ += total_flops;
00488 #endif
00489
00490 IsComputed_ = true;
00491 NumCompute_++;
00492 ComputeTime_ += Time_.ElapsedTime();
00493
00494 return(ierr);
00495
00496 }
00497
00498
00499
00500 int Ifpack_ILU::Solve(bool Trans, const Epetra_MultiVector& X,
00501 Epetra_MultiVector& Y) const
00502 {
00503
00504 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00505 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::ApplyInverse - Solve");
00506 #endif
00507
00508
00509 bool Upper = true;
00510 bool Lower = false;
00511 bool UnitDiagonal = true;
00512
00513 if (!Trans) {
00514
00515 IFPACK_CHK_ERR(L_->Solve(Lower, Trans, UnitDiagonal, X, Y));
00516
00517 IFPACK_CHK_ERR(Y.Multiply(1.0, *D_, Y, 0.0));
00518
00519 IFPACK_CHK_ERR(U_->Solve(Upper, Trans, UnitDiagonal, Y, Y));
00520 }
00521 else {
00522
00523 IFPACK_CHK_ERR(U_->Solve(Upper, Trans, UnitDiagonal, X, Y));
00524
00525 IFPACK_CHK_ERR(Y.Multiply(1.0, *D_, Y, 0.0));
00526 IFPACK_CHK_ERR(L_->Solve(Lower, Trans, UnitDiagonal, Y, Y));
00527 }
00528
00529
00530 return(0);
00531 }
00532
00533
00534
00535 int Ifpack_ILU::Multiply(bool Trans, const Epetra_MultiVector& X,
00536 Epetra_MultiVector& Y) const
00537 {
00538
00539 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00540 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::Multiply");
00541 #endif
00542
00543 if (!IsComputed())
00544 IFPACK_CHK_ERR(-3);
00545
00546 if (!Trans) {
00547 IFPACK_CHK_ERR(U_->Multiply(Trans, X, Y));
00548
00549 IFPACK_CHK_ERR(Y.Update(1.0, X, 1.0));
00550
00551 IFPACK_CHK_ERR(Y.ReciprocalMultiply(1.0, *D_, Y, 0.0));
00552 Epetra_MultiVector Y1temp(Y);
00553 IFPACK_CHK_ERR(L_->Multiply(Trans, Y1temp, Y));
00554
00555 IFPACK_CHK_ERR(Y.Update(1.0, Y1temp, 1.0));
00556 }
00557 else {
00558
00559 IFPACK_CHK_ERR(L_->Multiply(Trans, X, Y));
00560
00561 IFPACK_CHK_ERR(Y.Update(1.0, X, 1.0));
00562
00563 IFPACK_CHK_ERR(Y.ReciprocalMultiply(1.0, *D_, Y, 0.0));
00564 Epetra_MultiVector Y1temp(Y);
00565 IFPACK_CHK_ERR(U_->Multiply(Trans, Y1temp, Y));
00566
00567 IFPACK_CHK_ERR(Y.Update(1.0, Y1temp, 1.0));
00568 }
00569
00570 return(0);
00571 }
00572
00573
00574
00575 int Ifpack_ILU::ApplyInverse(const Epetra_MultiVector& X,
00576 Epetra_MultiVector& Y) const
00577 {
00578
00579 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00580 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::ApplyInverse");
00581 #endif
00582
00583 if (!IsComputed())
00584 IFPACK_CHK_ERR(-3);
00585
00586 if (X.NumVectors() != Y.NumVectors())
00587 IFPACK_CHK_ERR(-2);
00588
00589 Time_.ResetStartTime();
00590
00591
00592
00593 Teuchos::RefCountPtr< const Epetra_MultiVector > Xcopy;
00594 if (X.Pointers()[0] == Y.Pointers()[0])
00595 Xcopy = Teuchos::rcp( new Epetra_MultiVector(X) );
00596 else
00597 Xcopy = Teuchos::rcp( &X, false );
00598
00599 IFPACK_CHK_ERR(Solve(Ifpack_ILU::UseTranspose(), *Xcopy, Y));
00600
00601
00602 #ifdef IFPACK_FLOPCOUNTERS
00603 ApplyInverseFlops_ += X.NumVectors() * 4 *
00604 (L_->NumGlobalNonzeros() + U_->NumGlobalNonzeros());
00605 #endif
00606
00607 ++NumApplyInverse_;
00608 ApplyInverseTime_ += Time_.ElapsedTime();
00609
00610 return(0);
00611
00612 }
00613
00614
00615 double Ifpack_ILU::Condest(const Ifpack_CondestType CT,
00616 const int MaxIters, const double Tol,
00617 Epetra_RowMatrix* Matrix_in)
00618 {
00619
00620 #ifdef IFPACK_TEUCHOS_TIME_MONITOR
00621 TEUCHOS_FUNC_TIME_MONITOR("Ifpack_ILU::Condest");
00622 #endif
00623
00624 if (!IsComputed())
00625 return(-1.0);
00626
00627 Condest_ = Ifpack_Condest(*this, CT, MaxIters, Tol, Matrix_in);
00628
00629 return(Condest_);
00630 }
00631
00632
00633 std::ostream&
00634 Ifpack_ILU::Print(std::ostream& os) const
00635 {
00636 if (!Comm().MyPID()) {
00637 os << endl;
00638 os << "================================================================================" << endl;
00639 os << "Ifpack_ILU: " << Label() << endl << endl;
00640 os << "Level-of-fill = " << LevelOfFill() << endl;
00641 os << "Absolute threshold = " << AbsoluteThreshold() << endl;
00642 os << "Relative threshold = " << RelativeThreshold() << endl;
00643 os << "Relax value = " << RelaxValue() << endl;
00644 os << "Condition number estimate = " << Condest() << endl;
00645 os << "Global number of rows = " << A_->NumGlobalRows64() << endl;
00646 if (IsComputed_) {
00647 os << "Number of rows of L, D, U = " << L_->NumGlobalRows64() << endl;
00648 os << "Number of nonzeros of L + U = " << NumGlobalNonzeros64() << endl;
00649 os << "nonzeros / rows = "
00650 << 1.0 * NumGlobalNonzeros64() / U_->NumGlobalRows64() << endl;
00651 }
00652 os << endl;
00653 os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
00654 os << "----- ------- -------------- ------------ --------" << endl;
00655 os << "Initialize() " << std::setw(5) << NumInitialize()
00656 << " " << std::setw(15) << InitializeTime()
00657 << " 0.0 0.0" << endl;
00658 os << "Compute() " << std::setw(5) << NumCompute()
00659 << " " << std::setw(15) << ComputeTime()
00660 << " " << std::setw(15) << 1.0e-6 * ComputeFlops();
00661 if (ComputeTime() != 0.0)
00662 os << " " << std::setw(15) << 1.0e-6 * ComputeFlops() / ComputeTime() << endl;
00663 else
00664 os << " " << std::setw(15) << 0.0 << endl;
00665 os << "ApplyInverse() " << std::setw(5) << NumApplyInverse()
00666 << " " << std::setw(15) << ApplyInverseTime()
00667 << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops();
00668 if (ApplyInverseTime() != 0.0)
00669 os << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops() / ApplyInverseTime() << endl;
00670 else
00671 os << " " << std::setw(15) << 0.0 << endl;
00672 os << "================================================================================" << endl;
00673 os << endl;
00674 }
00675
00676 return(os);
00677 }