Ifpack_ICT.cpp

Go to the documentation of this file.
00001 /*@HEADER
00002 // ***********************************************************************
00003 // 
00004 //       Ifpack: Object-Oriented Algebraic Preconditioner Package
00005 //                 Copyright (2002) Sandia Corporation
00006 // 
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 // 
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //  
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //  
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 
00025 // 
00026 // ***********************************************************************
00027 //@HEADER
00028 */
00029 
00030 #include "Ifpack_ConfigDefs.h"
00031 #ifdef HAVE_IFPACK_TEUCHOS
00032 #include "Ifpack_Preconditioner.h"
00033 #include "Ifpack_ICT.h"
00034 #include "Ifpack_Condest.h"
00035 #include "Ifpack_Utils.h"
00036 #include "Ifpack_HashTable.h"
00037 #include "Epetra_SerialComm.h"
00038 #include "Epetra_Comm.h"
00039 #include "Epetra_Map.h"
00040 #include "Epetra_RowMatrix.h"
00041 #include "Epetra_CrsMatrix.h"
00042 #include "Epetra_Vector.h"
00043 #include "Epetra_MultiVector.h"
00044 #include "Epetra_Util.h"
00045 #include "Teuchos_ParameterList.hpp"
00046 #include "Teuchos_RefCountPtr.hpp"
00047 
00048 //==============================================================================
00049 // FIXME: allocate Comm_ and Time_ the first Initialize() call
00050 Ifpack_ICT::Ifpack_ICT(const Epetra_RowMatrix* A) :
00051   A_(*A),
00052   Comm_(A_.Comm()),
00053   Condest_(-1.0),
00054   Athresh_(0.0),
00055   Rthresh_(1.0),
00056   LevelOfFill_(1.0),
00057   DropTolerance_(0.0),
00058   Relax_(0.0),
00059   IsInitialized_(false),
00060   IsComputed_(false),
00061   UseTranspose_(false),
00062   NumMyRows_(0),
00063   NumInitialize_(0),
00064   NumCompute_(0),
00065   NumApplyInverse_(0),
00066   InitializeTime_(0.0),
00067   ComputeTime_(0.0),
00068   ApplyInverseTime_(0.0),
00069   ComputeFlops_(0.0),
00070   ApplyInverseFlops_(0.0),
00071   Time_(Comm()),
00072   GlobalNonzeros_(0)
00073 {
00074   // do nothing here
00075 }
00076 
00077 //==============================================================================
00078 Ifpack_ICT::~Ifpack_ICT()
00079 {
00080   Destroy();
00081 }
00082 
00083 //==============================================================================
00084 void Ifpack_ICT::Destroy()
00085 {
00086   IsInitialized_ = false;
00087   IsComputed_ = false;
00088 }
00089 
00090 //==========================================================================
00091 int Ifpack_ICT::SetParameters(Teuchos::ParameterList& List)
00092 {
00093 
00094   try
00095   {
00096     LevelOfFill_ = List.get("fact: ict level-of-fill",LevelOfFill_);
00097     Athresh_ = List.get("fact: absolute threshold", Athresh_);
00098     Rthresh_ = List.get("fact: relative threshold", Rthresh_);
00099     Relax_ = List.get("fact: relax value", Relax_);
00100     DropTolerance_ = List.get("fact: drop tolerance", DropTolerance_);
00101 
00102     // set label
00103     Label_ = "ICT (fill=" + Ifpack_toString(LevelOfFill())
00104       + ", athr=" + Ifpack_toString(AbsoluteThreshold()) 
00105       + ", rthr=" + Ifpack_toString(RelativeThreshold())
00106       + ", relax=" + Ifpack_toString(RelaxValue())
00107       + ", droptol=" + Ifpack_toString(DropTolerance())
00108       + ")";
00109 
00110     return(0);
00111   }
00112   catch (...)
00113   {
00114     cerr << "Caught an exception while parsing the parameter list" << endl;
00115     cerr << "This typically means that a parameter was set with the" << endl;
00116     cerr << "wrong type (for example, int instead of double). " << endl;
00117     cerr << "please check the documentation for the type required by each parameer." << endl;
00118     IFPACK_CHK_ERR(-1);
00119   }
00120 }
00121 
00122 //==========================================================================
00123 int Ifpack_ICT::Initialize()
00124 {
00125   // clean data if present
00126   Destroy();
00127 
00128   Time_.ResetStartTime();
00129 
00130   // matrix must be square. Check only on one processor
00131   if (Comm().NumProc() == 1 && Matrix().NumMyRows() != Matrix().NumMyCols())
00132     IFPACK_CHK_ERR(-2);
00133     
00134   NumMyRows_ = Matrix().NumMyRows();
00135 
00136   // nothing else to do here
00137   IsInitialized_ = true;
00138   ++NumInitialize_;
00139   InitializeTime_ += Time_.ElapsedTime();
00140 
00141   return(0);
00142 }
00143 
00144 //==========================================================================
00145 int Ifpack_ICT::Compute() 
00146 {
00147   if (!IsInitialized()) 
00148     IFPACK_CHK_ERR(Initialize());
00149 
00150   Time_.ResetStartTime();
00151   IsComputed_ = false;
00152 
00153   NumMyRows_ = A_.NumMyRows();
00154   int Length = A_.MaxNumEntries();
00155   vector<int>    RowIndices(Length);
00156   vector<double> RowValues(Length);
00157 
00158   bool distributed = (Comm().NumProc() > 1)?true:false;
00159 
00160   if (distributed)
00161   {
00162     SerialComm_ = Teuchos::rcp(new Epetra_SerialComm);
00163     SerialMap_ = Teuchos::rcp(new Epetra_Map(NumMyRows_, 0, *SerialComm_));
00164     assert (SerialComm_.get() != 0);
00165     assert (SerialMap_.get() != 0);
00166   }
00167   else
00168     SerialMap_ = Teuchos::rcp(const_cast<Epetra_Map*>(&A_.RowMatrixRowMap()), false);
00169 
00170   int RowNnz;
00171   double flops = 0.0;
00172 
00173   H_ = Teuchos::rcp(new Epetra_CrsMatrix(Copy,*SerialMap_,0));
00174   if (H_.get() == 0)
00175     IFPACK_CHK_ERR(-5); // memory allocation error
00176 
00177   // get A(0,0) element and insert it (after sqrt)
00178   IFPACK_CHK_ERR(A_.ExtractMyRowCopy(0,Length,RowNnz,
00179                                      &RowValues[0],&RowIndices[0]));
00180 
00181   // skip off-processor elements
00182   if (distributed)
00183   {
00184     int count = 0;
00185     for (int i = 0 ;i < RowNnz ; ++i) 
00186     {
00187       if (RowIndices[i] < NumMyRows_){
00188         RowIndices[count] = RowIndices[i];
00189         RowValues[count] = RowValues[i];
00190         ++count;
00191       }
00192       else
00193         continue;
00194     }
00195     RowNnz = count;
00196   }
00197 
00198   // modify diagonal
00199   double diag_val = 0.0;
00200   for (int i = 0 ;i < RowNnz ; ++i) {
00201     if (RowIndices[i] == 0) {
00202       double& v = RowValues[i];
00203       diag_val = AbsoluteThreshold() * EPETRA_SGN(v) +
00204         RelativeThreshold() * v;
00205       break;
00206     }
00207   }
00208 
00209   diag_val = sqrt(diag_val);
00210   int diag_idx = 0;
00211   EPETRA_CHK_ERR(H_->InsertGlobalValues(0,1,&diag_val, &diag_idx));
00212 
00213   int oldSize = RowNnz;
00214 
00215   // this is a bit magic
00216   int hash_size = 1024;
00217   while (hash_size < (int)(A_.MaxNumEntries() * LevelOfFill()))
00218     hash_size *= 2;
00219   Ifpack_HashTable Hash(hash_size - 1, 1);
00220 
00221   // start factorization for line 1
00222   for (int row_i = 1 ; row_i < NumMyRows_ ; ++row_i) {
00223 
00224     // get row `row_i' of the matrix
00225     IFPACK_CHK_ERR(A_.ExtractMyRowCopy(row_i,Length,RowNnz,
00226                                        &RowValues[0],&RowIndices[0]));
00227 
00228     // skip off-processor elements
00229     if (distributed)
00230     {
00231       int count = 0;
00232       for (int i = 0 ;i < RowNnz ; ++i) 
00233       {
00234         if (RowIndices[i] < NumMyRows_){
00235           RowIndices[count] = RowIndices[i];
00236           RowValues[count] = RowValues[i];
00237           ++count;
00238         }
00239         else
00240           continue;
00241       }
00242       RowNnz = count;
00243     }
00244 
00245     // number of nonzeros in this row are defined as the nonzeros
00246     // of the matrix, plus the level of fill 
00247     int LOF = (int)(LevelOfFill() * RowNnz);
00248     if (LOF == 0) LOF = 1;
00249 
00250     // convert line `row_i' into hash for fast access
00251     Hash.reset();
00252 
00253     double h_ii = 0.0;
00254     for (int i = 0 ; i < RowNnz ; ++i) {
00255       if (RowIndices[i] == row_i) {
00256         double& v = RowValues[i];
00257         h_ii = AbsoluteThreshold() * EPETRA_SGN(v) + RelativeThreshold() * v;
00258       }
00259       else if (RowIndices[i] < row_i)
00260       {
00261         Hash.set(RowIndices[i], RowValues[i], true);
00262       }
00263     }
00264       
00265     // form element (row_i, col_j)
00266     // I start from the first row that has a nonzero column
00267     // index in row_i.
00268     for (int col_j = RowIndices[0] ; col_j < row_i ; ++col_j) {
00269 
00270       short int flops = 0;
00271       double h_ij = 0.0, h_jj = 0.0;
00272       // note: get() returns 0.0 if col_j is not found
00273       h_ij = Hash.get(col_j);
00274 
00275       // get pointers to row `col_j'
00276       int* ColIndices;
00277       double* ColValues;
00278       int ColNnz;
00279       H_->ExtractGlobalRowView(col_j, ColNnz, ColValues, ColIndices);
00280 
00281       for (int k = 0 ; k < ColNnz ; ++k) {
00282         int col_k = ColIndices[k];
00283 
00284         if (col_k == col_j)
00285           h_jj = ColValues[k];
00286         else {
00287           double xxx = Hash.get(col_k);
00288           if (xxx != 0.0)
00289           {
00290             h_ij -= ColValues[k] * xxx;
00291             flops += 2;
00292           }
00293         }
00294       }
00295 
00296       h_ij /= h_jj;
00297 
00298       if (IFPACK_ABS(h_ij) > DropTolerance_)
00299       {
00300         Hash.set(col_j, h_ij);
00301       }
00302     
00303       // only approx
00304       ComputeFlops_ += 2.0 * flops + 1;
00305     }
00306 
00307     int size = Hash.getNumEntries();
00308 
00309     vector<double> AbsRow(size);
00310     int count = 0;
00311     
00312     // +1 because I use the extra position for diagonal in insert
00313     vector<int> keys(size + 1);
00314     vector<double> values(size + 1);
00315 
00316     Hash.arrayify(&keys[0], &values[0]);
00317 
00318     for (int i = 0 ; i < size ; ++i)
00319     {
00320       AbsRow[i] = IFPACK_ABS(values[i]);
00321     }
00322     count = size;
00323 
00324     double cutoff = 0.0;
00325     if (count > LOF) {
00326       nth_element(AbsRow.begin(), AbsRow.begin() + LOF, AbsRow.begin() + count, 
00327                   greater<double>());
00328       cutoff = AbsRow[LOF];
00329     }
00330 
00331     for (int i = 0 ; i < size ; ++i)
00332     {
00333       h_ii -= values[i] * values[i];
00334     }
00335 
00336     if (h_ii < 0.0) h_ii = 1e-12;;
00337 
00338     h_ii = sqrt(h_ii);
00339 
00340     // only approx, + 1 == sqrt
00341     ComputeFlops_ += 2 * size + 1;
00342 
00343     double DiscardedElements = 0.0;
00344 
00345     count = 0;
00346     for (int i = 0 ; i < size ; ++i)    
00347     { 
00348       if (IFPACK_ABS(values[i]) > cutoff)
00349       {
00350         values[count] = values[i];
00351         keys[count] = keys[i];
00352         ++count;
00353       }
00354       else  
00355         DiscardedElements += values[i];
00356     }
00357 
00358     if (RelaxValue() != 0.0) {
00359       DiscardedElements *= RelaxValue();
00360       h_ii += DiscardedElements;
00361     }
00362 
00363     values[count] = h_ii;
00364     keys[count] = row_i;
00365     ++count;
00366 
00367     H_->InsertGlobalValues(row_i, count, &(values[0]), (int*)&(keys[0]));
00368 
00369     oldSize = size;
00370   }
00371 
00372   IFPACK_CHK_ERR(H_->FillComplete());
00373 
00374 #if 0
00375   // to check the complete factorization
00376   Epetra_Vector LHS(Matrix().RowMatrixRowMap());
00377   Epetra_Vector RHS1(Matrix().RowMatrixRowMap());
00378   Epetra_Vector RHS2(Matrix().RowMatrixRowMap());
00379   Epetra_Vector RHS3(Matrix().RowMatrixRowMap());
00380   LHS.Random();
00381 
00382   Matrix().Multiply(false,LHS,RHS1);
00383   H_->Multiply(true,LHS,RHS2);
00384   H_->Multiply(false,RHS2,RHS3);
00385 
00386   RHS1.Update(-1.0, RHS3, 1.0);
00387   cout << endl;
00388   cout << RHS1;
00389 #endif
00390   int MyNonzeros = H_->NumGlobalNonzeros();
00391   Comm().SumAll(&MyNonzeros, &GlobalNonzeros_, 1);
00392 
00393   IsComputed_ = true;
00394   double TotalFlops; // sum across all the processors
00395   A_.Comm().SumAll(&flops, &TotalFlops, 1);
00396   ComputeFlops_ += TotalFlops;
00397   ++NumCompute_;
00398   ComputeTime_ += Time_.ElapsedTime();
00399 
00400   return(0);
00401 
00402 }
00403 
00404 //=============================================================================
00405 int Ifpack_ICT::ApplyInverse(const Epetra_MultiVector& X, 
00406            Epetra_MultiVector& Y) const
00407 {
00408 
00409   if (!IsComputed())
00410     IFPACK_CHK_ERR(-3); // compute preconditioner first
00411 
00412   if (X.NumVectors() != Y.NumVectors()) 
00413     IFPACK_CHK_ERR(-2); // Return error: X and Y not the same size
00414 
00415   Time_.ResetStartTime();
00416 
00417   // AztecOO gives X and Y pointing to the same memory location,
00418   // need to create an auxiliary vector, Xcopy
00419   const Epetra_MultiVector* Xcopy;
00420   if (X.Pointers()[0] == Y.Pointers()[0])
00421     Xcopy = new Epetra_MultiVector(X);
00422   else
00423     Xcopy = &X;
00424 
00425   // NOTE: H_ is based on SerialMap_, while Xcopy is based
00426   // on A.Map()... which are in general different. However, Solve()
00427   // does not seem to care... which is fine with me.
00428   //
00429   EPETRA_CHK_ERR(H_->Solve(false,false,false,*Xcopy,Y));
00430   EPETRA_CHK_ERR(H_->Solve(false,true,false,Y,Y));
00431 
00432   if (Xcopy != &X)
00433     delete Xcopy;
00434 
00435   // these are global flop count
00436   ApplyInverseFlops_ += 4.0 * GlobalNonzeros_;
00437 
00438   ++NumApplyInverse_;
00439   ApplyInverseTime_ += Time_.ElapsedTime();
00440 
00441   return(0);
00442 }
00443 //=============================================================================
00444 // This function finds X such that LDU Y = X or U(trans) D L(trans) Y = X for multiple RHS
00445 int Ifpack_ICT::Apply(const Epetra_MultiVector& X, 
00446           Epetra_MultiVector& Y) const 
00447 {
00448 
00449   IFPACK_CHK_ERR(-98);
00450 }
00451 
00452 //=============================================================================
00453 double Ifpack_ICT::Condest(const Ifpack_CondestType CT, 
00454                             const int MaxIters, const double Tol,
00455           Epetra_RowMatrix* Matrix)
00456 {
00457   if (!IsComputed()) // cannot compute right now
00458     return(-1.0);
00459 
00460   // NOTE: this is computing the *local* condest
00461   if (Condest_ == -1.0)
00462     Condest_ = Ifpack_Condest(*this, CT, MaxIters, Tol, Matrix);
00463 
00464   return(Condest_);
00465 }
00466 
00467 //=============================================================================
00468 std::ostream&
00469 Ifpack_ICT::Print(std::ostream& os) const
00470 {
00471   if (!Comm().MyPID()) {
00472     os << endl;
00473     os << "================================================================================" << endl;
00474     os << "Ifpack_ICT: " << Label() << endl << endl;
00475     os << "Level-of-fill      = " << LevelOfFill() << endl;
00476     os << "Absolute threshold = " << AbsoluteThreshold() << endl;
00477     os << "Relative threshold = " << RelativeThreshold() << endl;
00478     os << "Relax value        = " << RelaxValue() << endl;
00479     os << "Condition number estimate = " << Condest() << endl;
00480     os << "Global number of rows            = " << Matrix().NumGlobalRows() << endl;
00481     if (IsComputed_) {
00482       os << "Number of nonzeros of H         = " << H_->NumGlobalNonzeros() << endl;
00483       os << "nonzeros / rows                 = " 
00484          << 1.0 * H_->NumGlobalNonzeros() / H_->NumGlobalRows() << endl;
00485     }
00486     os << endl;
00487     os << "Phase           # calls   Total Time (s)       Total MFlops     MFlops/s" << endl;
00488     os << "-----           -------   --------------       ------------     --------" << endl;
00489     os << "Initialize()    "   << std::setw(5) << NumInitialize() 
00490        << "  " << std::setw(15) << InitializeTime() 
00491        << "               0.0            0.0" << endl;
00492     os << "Compute()       "   << std::setw(5) << NumCompute() 
00493        << "  " << std::setw(15) << ComputeTime()
00494        << "  " << std::setw(15) << 1.0e-6 * ComputeFlops();
00495     if (ComputeTime() != 0.0) 
00496       os << "  " << std::setw(15) << 1.0e-6 * ComputeFlops() / ComputeTime() << endl;
00497     else
00498       os << "  " << std::setw(15) << 0.0 << endl;
00499     os << "ApplyInverse()  "   << std::setw(5) << NumApplyInverse() 
00500        << "  " << std::setw(15) << ApplyInverseTime()
00501        << "  " << std::setw(15) << 1.0e-6 * ApplyInverseFlops();
00502     if (ApplyInverseTime() != 0.0)
00503       os << "  " << std::setw(15) << 1.0e-6 * ApplyInverseFlops() / ApplyInverseTime() << endl;
00504     else
00505       os << "  " << std::setw(15) << 0.0 << endl;
00506     os << "================================================================================" << endl;
00507     os << endl;
00508   }
00509 
00510   
00511   return(os);
00512 }
00513 #endif

Generated on Thu Sep 18 12:37:21 2008 for Ifpack Package Browser (Single Doxygen Collection) by doxygen 1.3.9.1