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_IC.h"
00046 #include "Ifpack_IC_Utils.h"
00047 #include "Ifpack_Condest.h"
00048 #include "Epetra_Comm.h"
00049 #include "Epetra_Map.h"
00050 #include "Epetra_RowMatrix.h"
00051 #include "Epetra_CrsMatrix.h"
00052 #include "Epetra_Vector.h"
00053 #include "Epetra_MultiVector.h"
00054 #include "Epetra_Util.h"
00055 #include "Teuchos_ParameterList.hpp"
00056 #include "Teuchos_RefCountPtr.hpp"
00057 using Teuchos::RefCountPtr;
00058 using Teuchos::rcp;
00059
00060
00061 Ifpack_IC::Ifpack_IC(Epetra_RowMatrix* A) :
00062 A_(rcp(A,false)),
00063 Comm_(A->Comm()),
00064 UseTranspose_(false),
00065 Condest_(-1.0),
00066 Athresh_(0.0),
00067 Rthresh_(1.0),
00068 Droptol_(0.0),
00069 Lfil_(1.0),
00070 Aict_(0),
00071 Lict_(0),
00072 Ldiag_(0),
00073 IsInitialized_(false),
00074 IsComputed_(false),
00075 NumInitialize_(0),
00076 NumCompute_(0),
00077 NumApplyInverse_(0),
00078 InitializeTime_(0.0),
00079 ComputeTime_(0),
00080 ApplyInverseTime_(0),
00081 Time_(Comm()),
00082 ComputeFlops_(0.0),
00083 ApplyInverseFlops_(0.0)
00084 {
00085 Teuchos::ParameterList List;
00086 SetParameters(List);
00087
00088 }
00089
00090 Ifpack_IC::~Ifpack_IC()
00091 {
00092 if (Lict_ != 0) {
00093 Ifpack_AIJMatrix * Lict = (Ifpack_AIJMatrix *) Lict_;
00094 delete [] Lict->ptr;
00095 delete [] Lict->col;
00096 delete [] Lict->val;
00097 delete Lict;
00098 }
00099 if (Aict_ != 0) {
00100 Ifpack_AIJMatrix * Aict = (Ifpack_AIJMatrix *) Aict_;
00101 delete Aict;
00102 }
00103 if (Ldiag_ != 0) delete [] Ldiag_;
00104
00105 IsInitialized_ = false;
00106 IsComputed_ = false;
00107 }
00108
00109
00110 int Ifpack_IC::SetParameters(Teuchos::ParameterList& List)
00111 {
00112
00113
00114 Lfil_ = List.get("fact: ict level-of-fill", Lfil_);
00115 Athresh_ = List.get("fact: absolute threshold", Athresh_);
00116 Rthresh_ = List.get("fact: relative threshold", Rthresh_);
00117 Droptol_ = List.get("fact: drop tolerance", Droptol_);
00118
00119
00120 sprintf(Label_, "IFPACK IC (fill=%f, drop=%f)",
00121 Lfil_, Droptol_);
00122 return(0);
00123 }
00124
00125
00126 int Ifpack_IC::Initialize()
00127 {
00128
00129 IsInitialized_ = false;
00130
00131
00132
00133 ++NumInitialize_;
00134 IsInitialized_ = true;
00135 return(0);
00136 }
00137
00138
00139 int Ifpack_IC::ComputeSetup()
00140 {
00141
00142 U_ = rcp(new Epetra_CrsMatrix(Copy, Matrix().RowMatrixRowMap(),
00143 Matrix().RowMatrixRowMap(), 0));
00144 D_ = rcp(new Epetra_Vector(Matrix().RowMatrixRowMap()));
00145
00146 if (U_.get() == 0 || D_.get() == 0)
00147 IFPACK_CHK_ERR(-5);
00148
00149 int ierr = 0;
00150 int i, j;
00151 int NumIn, NumU;
00152 bool DiagFound;
00153 int NumNonzeroDiags = 0;
00154
00155
00156 int MaxNumEntries = Matrix().MaxNumEntries();
00157
00158 std::vector<int> InI(MaxNumEntries);
00159 std::vector<int> UI(MaxNumEntries);
00160 std::vector<double> InV(MaxNumEntries);
00161 std::vector<double> UV(MaxNumEntries);
00162
00163 double *DV;
00164 ierr = D_->ExtractView(&DV);
00165
00166
00167
00168 int NumRows = Matrix().NumMyRows();
00169
00170 for (i=0; i< NumRows; i++) {
00171
00172 Matrix().ExtractMyRowCopy(i, MaxNumEntries, NumIn, &InV[0], &InI[0]);
00173
00174
00175 NumU = 0;
00176 DiagFound = false;
00177
00178 for (j=0; j< NumIn; j++) {
00179 int k = InI[j];
00180
00181 if (k==i) {
00182 DiagFound = true;
00183
00184 DV[i] += Rthresh_ * InV[j] + EPETRA_SGN(InV[j]) * Athresh_;
00185 }
00186
00187 else if (k < 0) return(-1);
00188 else if (i<k && k<NumRows) {
00189 UI[NumU] = k;
00190 UV[NumU] = InV[j];
00191 NumU++;
00192 }
00193 }
00194
00195
00196
00197 if (DiagFound) NumNonzeroDiags++;
00198 if (NumU) U_->InsertMyValues(i, NumU, &UV[0], &UI[0]);
00199
00200 }
00201
00202 U_->FillComplete(Matrix().OperatorDomainMap(),
00203 Matrix().OperatorRangeMap());
00204
00205 int ierr1 = 0;
00206 if (NumNonzeroDiags<U_->NumMyRows()) ierr1 = 1;
00207 Matrix().Comm().MaxAll(&ierr1, &ierr, 1);
00208 IFPACK_CHK_ERR(ierr);
00209
00210 return(0);
00211 }
00212
00213
00214 int Ifpack_IC::Compute() {
00215
00216 if (!IsInitialized())
00217 IFPACK_CHK_ERR(Initialize());
00218
00219 Time_.ResetStartTime();
00220 IsComputed_ = false;
00221
00222
00223 IFPACK_CHK_ERR(ComputeSetup());
00224
00225 int i;
00226
00227 int m, n, nz, Nrhs, ldrhs, ldlhs;
00228 int * ptr=0, * ind;
00229 double * val, * rhs, * lhs;
00230
00231 int ierr = Epetra_Util_ExtractHbData(U_.get(), 0, 0, m, n, nz, ptr, ind,
00232 val, Nrhs, rhs, ldrhs, lhs, ldlhs);
00233 if (ierr < 0)
00234 IFPACK_CHK_ERR(ierr);
00235
00236 Ifpack_AIJMatrix * Aict;
00237 if (Aict_==0) {
00238 Aict = new Ifpack_AIJMatrix;
00239 Aict_ = (void *) Aict;
00240 }
00241 else Aict = (Ifpack_AIJMatrix *) Aict_;
00242 Ifpack_AIJMatrix * Lict;
00243 if (Lict_==0) {
00244 Lict = new Ifpack_AIJMatrix;
00245 Lict_ = (void *) Lict;
00246 }
00247 else {
00248 Lict = (Ifpack_AIJMatrix *) Lict_;
00249 Ifpack_AIJMatrix_dealloc( Lict );
00250 }
00251 if (Ldiag_ != 0) delete [] Ldiag_;
00252 Aict->val = val;
00253 Aict->col = ind;
00254 Aict->ptr = ptr;
00255 double *DV;
00256 EPETRA_CHK_ERR(D_->ExtractView(&DV));
00257
00258
00259
00260
00261 int lfil = (Lfil_ * nz)/n + .5;
00262
00263 crout_ict(m, Aict, DV, Droptol_, lfil, Lict, &Ldiag_);
00264
00265
00266 delete [] ptr;
00267
00268
00269 U_ = rcp(new Epetra_CrsMatrix(View, A_->RowMatrixRowMap(), A_->RowMatrixRowMap(),0));
00270 D_ = rcp(new Epetra_Vector(View, A_->RowMatrixRowMap(), Ldiag_));
00271
00272 ptr = Lict->ptr;
00273 ind = Lict->col;
00274 val = Lict->val;
00275
00276 for (i=0; i< m; i++) {
00277 int NumEntries = ptr[i+1]-ptr[i];
00278 int * Indices = ind+ptr[i];
00279 double * Values = val+ptr[i];
00280 U_->InsertMyValues(i, NumEntries, Values, Indices);
00281 }
00282
00283 U_->FillComplete(A_->OperatorDomainMap(), A_->OperatorRangeMap());
00284 D_->Reciprocal(*D_);
00285
00286 #ifdef IFPACK_FLOPCOUNTERS
00287 double current_flops = 2 * nz;
00288 double total_flops = 0;
00289
00290 A_->Comm().SumAll(¤t_flops, &total_flops, 1);
00291
00292 ComputeFlops_ += total_flops;
00293
00294 ComputeFlops_ += (double) U_->NumGlobalNonzeros();
00295 ComputeFlops_ += (double) D_->GlobalLength();
00296 #endif
00297 ++NumCompute_;
00298 ComputeTime_ += Time_.ElapsedTime();
00299
00300
00301 IsComputed_ = true;
00302
00303 return(0);
00304
00305 }
00306
00307
00308
00309 int Ifpack_IC::ApplyInverse(const Epetra_MultiVector& X,
00310 Epetra_MultiVector& Y) const
00311 {
00312
00313 if (!IsComputed())
00314 IFPACK_CHK_ERR(-3);
00315
00316 if (X.NumVectors() != Y.NumVectors())
00317 IFPACK_CHK_ERR(-2);
00318
00319 Time_.ResetStartTime();
00320
00321 bool Upper = true;
00322 bool UnitDiagonal = true;
00323
00324
00325
00326 RefCountPtr< const Epetra_MultiVector > Xcopy;
00327 if (X.Pointers()[0] == Y.Pointers()[0])
00328 Xcopy = rcp( new Epetra_MultiVector(X) );
00329 else
00330 Xcopy = rcp( &X, false );
00331
00332 U_->Solve(Upper, true, UnitDiagonal, *Xcopy, Y);
00333 Y.Multiply(1.0, *D_, Y, 0.0);
00334 U_->Solve(Upper, false, UnitDiagonal, Y, Y);
00335
00336 #ifdef IFPACK_FLOPCOUNTERS
00337 ApplyInverseFlops_ += 4.0 * U_->NumGlobalNonzeros();
00338 ApplyInverseFlops_ += D_->GlobalLength();
00339 #endif
00340
00341 ++NumApplyInverse_;
00342 ApplyInverseTime_ += Time_.ElapsedTime();
00343
00344 return(0);
00345
00346 }
00347
00348
00349
00350 int Ifpack_IC::Apply(const Epetra_MultiVector& X,
00351 Epetra_MultiVector& Y) const
00352 {
00353
00354 if (X.NumVectors() != Y.NumVectors())
00355 IFPACK_CHK_ERR(-2);
00356
00357 Epetra_MultiVector * X1 = (Epetra_MultiVector *) &X;
00358 Epetra_MultiVector * Y1 = (Epetra_MultiVector *) &Y;
00359
00360 U_->Multiply(false, *X1, *Y1);
00361 Y1->Update(1.0, *X1, 1.0);
00362 Y1->ReciprocalMultiply(1.0, *D_, *Y1, 0.0);
00363 Epetra_MultiVector Y1temp(*Y1);
00364 U_->Multiply(true, Y1temp, *Y1);
00365 Y1->Update(1.0, Y1temp, 1.0);
00366 return(0);
00367 }
00368
00369
00370 double Ifpack_IC::Condest(const Ifpack_CondestType CT,
00371 const int MaxIters, const double Tol,
00372 Epetra_RowMatrix* Matrix_in)
00373 {
00374 if (!IsComputed())
00375 return(-1.0);
00376
00377 if (Condest_ == -1.0)
00378 Condest_ = Ifpack_Condest(*this, CT, MaxIters, Tol, Matrix_in);
00379
00380 return(Condest_);
00381 }
00382
00383
00384 std::ostream&
00385 Ifpack_IC::Print(std::ostream& os) const
00386 {
00387 if (!Comm().MyPID()) {
00388 os << endl;
00389 os << "================================================================================" << endl;
00390 os << "Ifpack_IC: " << Label() << endl << endl;
00391 os << "Level-of-fill = " << LevelOfFill() << endl;
00392 os << "Absolute threshold = " << AbsoluteThreshold() << endl;
00393 os << "Relative threshold = " << RelativeThreshold() << endl;
00394 os << "Drop tolerance = " << DropTolerance() << endl;
00395 os << "Condition number estimate = " << Condest() << endl;
00396 os << "Global number of rows = " << A_->NumGlobalRows64() << endl;
00397 if (IsComputed_) {
00398 os << "Number of nonzeros of H = " << U_->NumGlobalNonzeros64() << endl;
00399 os << "nonzeros / rows = "
00400 << 1.0 * U_->NumGlobalNonzeros64() / U_->NumGlobalRows64() << endl;
00401 }
00402 os << endl;
00403 os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << endl;
00404 os << "----- ------- -------------- ------------ --------" << endl;
00405 os << "Initialize() " << std::setw(5) << NumInitialize()
00406 << " " << std::setw(15) << InitializeTime()
00407 << " 0.0 0.0" << endl;
00408 os << "Compute() " << std::setw(5) << NumCompute()
00409 << " " << std::setw(15) << ComputeTime()
00410 << " " << std::setw(15) << 1.0e-6 * ComputeFlops();
00411 if (ComputeTime() != 0.0)
00412 os << " " << std::setw(15) << 1.0e-6 * ComputeFlops() / ComputeTime() << endl;
00413 else
00414 os << " " << std::setw(15) << 0.0 << endl;
00415 os << "ApplyInverse() " << std::setw(5) << NumApplyInverse()
00416 << " " << std::setw(15) << ApplyInverseTime()
00417 << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops();
00418 if (ApplyInverseTime() != 0.0)
00419 os << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops() / ApplyInverseTime() << endl;
00420 else
00421 os << " " << std::setw(15) << 0.0 << endl;
00422 os << "================================================================================" << endl;
00423 os << endl;
00424 }
00425
00426
00427 return(os);
00428 }