Ifpack Package Browser (Single Doxygen Collection) Development
Ifpack_CrsRick.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_CrsRick.h"
00031 #include "Epetra_Comm.h"
00032 #include "Epetra_Map.h"
00033 #include "Epetra_CrsGraph.h"
00034 #include "Epetra_CrsMatrix.h"
00035 #include "Epetra_Vector.h"
00036 #include "Epetra_MultiVector.h"
00037 
00038 #include <Teuchos_ParameterList.hpp>
00039 #include <ifp_parameters.h>
00040 
00041 //==============================================================================
00042 Ifpack_CrsRick::Ifpack_CrsRick(const Epetra_CrsMatrix &A, const Ifpack_IlukGraph & Graph) 
00043   : A_(A),
00044     Graph_(Graph),
00045     UseTranspose_(false),
00046     Allocated_(false),
00047     ValuesInitialized_(false),
00048     Factored_(false),
00049     RelaxValue_(0.0),
00050     Condest_(-1.0),
00051     Athresh_(0.0),
00052     Rthresh_(1.0),
00053     OverlapX_(0),
00054     OverlapY_(0),
00055     OverlapMode_(Zero)
00056 {
00057   int ierr = Allocate();
00058 }
00059 
00060 //==============================================================================
00061 Ifpack_CrsRick::Ifpack_CrsRick(const Ifpack_CrsRick & FactoredMatrix) 
00062   : A_(FactoredMatrix.A_),
00063     Graph_(FactoredMatrix.Graph_),
00064     UseTranspose_(FactoredMatrix.UseTranspose_),
00065     Allocated_(FactoredMatrix.Allocated_),
00066     ValuesInitialized_(FactoredMatrix.ValuesInitialized_),
00067     Factored_(FactoredMatrix.Factored_),
00068     RelaxValue_(FactoredMatrix.RelaxValue_),
00069     Condest_(FactoredMatrix.Condest_),
00070     Athresh_(FactoredMatrix.Athresh_),
00071     Rthresh_(FactoredMatrix.Rthresh_),
00072     OverlapX_(0),
00073     OverlapY_(0),
00074     OverlapMode_(FactoredMatrix.OverlapMode_)
00075 {
00076   U_ = new Epetra_CrsMatrix(FactoredMatrix.U());
00077   D_ = new Epetra_Vector(Graph_.L_Graph().RowMap());
00078   
00079 }
00080 
00081 //==============================================================================
00082 int Ifpack_CrsRick::Allocate() {
00083 
00084   // Allocate Epetra_CrsMatrix using ILUK graphs
00085   U_ = new Epetra_CrsMatrix(Copy, Graph_.U_Graph());
00086   D_ = new Epetra_Vector(Graph_.L_Graph().RowMap());
00087   
00088   
00089     SetAllocated(true);
00090     return(0);
00091 }
00092 //==============================================================================
00093 Ifpack_CrsRick::~Ifpack_CrsRick(){
00094 
00095 
00096   delete U_;
00097   delete D_; // Diagonal is stored separately.  We store the inverse.
00098 
00099   if (OverlapX_!=0) delete OverlapX_;
00100   if (OverlapY_!=0) delete OverlapY_;
00101 
00102   ValuesInitialized_ = false;
00103   Factored_ = false;
00104   Allocated_ = false;
00105 }
00106 
00107 //==========================================================================
00108 int Ifpack_CrsRick::SetParameters(const Teuchos::ParameterList& parameterlist,
00109           bool cerr_warning_if_unused)
00110 {
00111   Ifpack::param_struct params;
00112   params.double_params[Ifpack::relax_value] = RelaxValue_;
00113   params.double_params[Ifpack::absolute_threshold] = Athresh_;
00114   params.double_params[Ifpack::relative_threshold] = Rthresh_;
00115   params.overlap_mode = OverlapMode_;
00116 
00117   Ifpack::set_parameters(parameterlist, params, cerr_warning_if_unused);
00118 
00119   RelaxValue_ = params.double_params[Ifpack::relax_value];
00120   Athresh_ = params.double_params[Ifpack::absolute_threshold];
00121   Rthresh_ = params.double_params[Ifpack::relative_threshold];
00122   OverlapMode_ = params.overlap_mode;
00123 
00124   return(0);
00125 }
00126 
00127 //==========================================================================
00128 int Ifpack_CrsRick::InitValues() {
00129 
00130   // if (!Allocated()) return(-1); // This test is not needed at this time.  All constructors allocate.
00131 
00132   int ierr = 0;
00133   int i, j;
00134   int * InI=0, * LI=0, * UI = 0;
00135   double * InV=0, * LV=0, * UV = 0;
00136   int NumIn, NumL, NumU;
00137   bool DiagFound;
00138   int NumNonzeroDiags = 0;
00139 
00140   Epetra_CrsMatrix * OverlapA = (Epetra_CrsMatrix *) &A_;
00141 
00142   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal()) {
00143   
00144   OverlapA = new Epetra_CrsMatrix(Copy, *Graph_.OverlapGraph());
00145   OverlapA->Import(A_, *Graph_.OverlapImporter(), Insert);
00146   }
00147 
00148   // Get Maximun Row length
00149   int MaxNumEntries = OverlapA->MaxNumEntries();
00150 
00151   InI = new int[MaxNumEntries]; // Allocate temp space
00152   UI = new int[MaxNumEntries];
00153   InV = new double[MaxNumEntries];
00154   UV = new double[MaxNumEntries];
00155 
00156   double *DV;
00157   ierr = D_->ExtractView(&DV); // Get view of diagonal
00158     
00159 
00160   // First we copy the user's matrix into diagonal vector and U, regardless of fill level
00161 
00162   for (i=0; i< NumMyRows(); i++) {
00163 
00164     OverlapA->ExtractMyRowCopy(i, MaxNumEntries, NumIn, InV, InI); // Get Values and Indices
00165     
00166     // Split into L and U (we don't assume that indices are ordered).
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_; // Store perturbed diagonal in Epetra_Vector D_
00178       }
00179 
00180       else if (k < 0) return(-1); // Out of range
00181       else if (k<NumMyRows()) {
00182   UI[NumU] = k;
00183   UV[NumU] = InV[j];
00184   NumU++;
00185       }
00186     }
00187     
00188     // Check in things for this row of L and U
00189 
00190     if (DiagFound) NumNonzeroDiags++;
00191     if (NumU) U_->ReplaceMyValues(i, NumU, UV, UI);
00192     
00193   }
00194 
00195   delete [] UI;
00196   delete [] UV;
00197   delete [] InI;
00198   delete [] InV;
00199 
00200   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal()) delete OverlapA;
00201 
00202 
00203   U_->FillComplete();
00204 
00205   // At this point L and U have the values of A in the structure of L and U, and diagonal vector D
00206 
00207   SetValuesInitialized(true);
00208   SetFactored(false);
00209 
00210   int TotalNonzeroDiags = 0;
00211   Graph_.U_Graph().RowMap().Comm().SumAll(&NumNonzeroDiags, &TotalNonzeroDiags, 1);
00212   if (Graph_.LevelOverlap()==0 &&
00213       ((TotalNonzeroDiags!=NumGlobalRows()) || 
00214        (TotalNonzeroDiags!=NumGlobalDiagonals()))) ierr = 1;
00215   if (NumNonzeroDiags != NumMyDiagonals()) ierr = 1; // Diagonals are not right, warn user
00216 
00217   return(ierr);
00218 }
00219 
00220 //==========================================================================
00221 int Ifpack_CrsRick::Factor() {
00222 
00223   // if (!Allocated()) return(-1); // This test is not needed at this time.  All constructors allocate.
00224   if (!ValuesInitialized()) return(-2); // Must have values initialized.
00225   if (Factored()) return(-3); // Can't have already computed factors.
00226 
00227   SetValuesInitialized(false);
00228 
00229   // MinMachNum should be officially defined, for now pick something a little 
00230   // bigger than IEEE underflow value
00231 
00232   double MinDiagonalValue = Epetra_MinDouble;
00233   double MaxDiagonalValue = 1.0/MinDiagonalValue;
00234 
00235   int ierr = 0;
00236   int i, j, k;
00237   int * UI = 0;
00238   double * UV = 0;
00239   int NumIn, NumU;
00240 
00241   // Get Maximun Row length
00242   int MaxNumEntries = U_->MaxNumEntries() + 1;
00243 
00244   int * InI = new int[MaxNumEntries]; // Allocate temp space
00245   double * InV = new double[MaxNumEntries];
00246   int * colflag = new int[NumMyCols()];
00247 
00248   double *DV;
00249   ierr = D_->ExtractView(&DV); // Get view of diagonal
00250 
00251   int current_madds = 0; // We will count multiply-add as they happen
00252 
00253   // Now start the factorization.
00254 
00255   // Need some integer workspace and pointers
00256   int NumUU; 
00257   int * UUI;
00258   double * UUV;
00259   for (j=0; j<NumMyCols(); j++) colflag[j] = - 1;
00260 
00261   for(i=0; i<NumMyRows(); i++) {
00262 
00263  // Fill InV, InI with current row of L, D and U combined
00264 
00265     NumIn = MaxNumEntries;
00266     IFPACK_CHK_ERR(L_->ExtractMyRowCopy(i, NumIn, NumL, InV, InI)==0);
00267     LV = InV;
00268     LI = InI;
00269 
00270     InV[NumL] = DV[i]; // Put in diagonal
00271     InI[NumL] = i;
00272     
00273     IFPACK_CHK_ERR(U_->ExtractMyRowCopy(i, NumIn-NumL-1, NumU, InV+NumL+1, InI+NumL+1));
00274     NumIn = NumL+NumU+1;
00275     UV = InV+NumL+1;
00276     UI = InI+NumL+1;
00277 
00278     // Set column flags
00279     for (j=0; j<NumIn; j++) colflag[InI[j]] = j;
00280 
00281     double diagmod = 0.0; // Off-diagonal accumulator
00282 
00283     for (int jj=0; jj<NumL; jj++) {
00284       j = InI[jj];
00285       double multiplier = InV[jj]; // current_mults++;
00286 
00287       InV[jj] *= DV[j];
00288       
00289       IFPACK_CHK_ERR(U_->ExtractMyRowView(j, NumUU, UUV, UUI)); // View of row above
00290 
00291       if (RelaxValue_==0.0) {
00292   for (k=0; k<NumUU; k++) {
00293     int kk = colflag[UUI[k]];
00294     if (kk>-1) {
00295       InV[kk] -= multiplier*UUV[k];
00296       current_madds++;
00297     }
00298   }
00299       }
00300       else {
00301   for (k=0; k<NumUU; k++) {
00302     int kk = colflag[UUI[k]];
00303     if (kk>-1) InV[kk] -= multiplier*UUV[k];
00304     else diagmod -= multiplier*UUV[k];
00305     current_madds++;
00306   }
00307       }
00308      }
00309     if (NumL) 
00310       IFPACK_CHK_ERR(L_->ReplaceMyValues(i, NumL, LV, LI));  // Replace current row of L
00311 
00312     DV[i] = InV[NumL]; // Extract Diagonal value
00313 
00314     if (RelaxValue_!=0.0) {
00315       DV[i] += RelaxValue_*diagmod; // Add off diagonal modifications
00316       // current_madds++;
00317     }
00318 
00319     if (fabs(DV[i]) > MaxDiagonalValue) {
00320       if (DV[i] < 0) DV[i] = - MinDiagonalValue;
00321       else DV[i] = MinDiagonalValue;
00322     }
00323     else
00324       DV[i] = 1.0/DV[i]; // Invert diagonal value
00325 
00326     for (j=0; j<NumU; j++) UV[j] *= DV[i]; // Scale U by inverse of diagonal
00327 
00328     if (NumU) 
00329       IFPACK_CHK_ERR(U_->ReplaceMyValues(i, NumU, UV, UI));  // Replace current row of L and U
00330 
00331 
00332     // Reset column flags
00333     for (j=0; j<NumIn; j++) colflag[InI[j]] = -1;
00334   }
00335 
00336   
00337   // Add up flops
00338  
00339   double current_flops = 2 * current_madds;
00340   double total_flops = 0;
00341     
00342   Graph_.L_Graph().RowMap().Comm().SumAll(&current_flops, &total_flops, 1); // Get total madds across all PEs
00343 
00344   // Now count the rest
00345   total_flops += (double) L_->NumGlobalNonzeros(); // Accounts for multiplier above
00346   total_flops += (double) D_->GlobalLength(); // Accounts for reciprocal of diagonal
00347   if (RelaxValue_!=0.0) total_flops += 2 * (double)D_->GlobalLength(); // Accounts for relax update of diag
00348 
00349   UpdateFlops(total_flops); // Update flop count
00350 
00351   delete [] InI;
00352   delete [] InV;
00353   delete [] colflag;
00354   
00355   SetFactored(true);
00356 
00357   return(ierr);
00358 
00359 }
00360 
00361 //=============================================================================
00362 int Ifpack_CrsRick::Solve(bool Trans, const Epetra_Vector& X, 
00363         Epetra_Vector& Y) const {
00364 //
00365 // This function finds Y such that LDU Y = X or U(trans) D L(trans) Y = X for a single RHS
00366 //
00367 
00368   bool Upper = true;
00369   bool Lower = false;
00370   bool UnitDiagonal = true;
00371 
00372   Epetra_Vector * X1 = (Epetra_Vector *) &X;
00373   Epetra_Vector * Y1 = (Epetra_Vector *) &Y;
00374 
00375   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal()) {
00376     if (OverlapX_==0) { // Need to allocate space for overlap X and Y
00377       OverlapX_ = new Epetra_Vector(Graph_.OverlapGraph()->RowMap());
00378       OverlapY_ = new Epetra_Vector(Graph_.OverlapGraph()->RowMap());
00379     }
00380     OverlapX_->Import(X,*Graph_.OverlapImporter(), Insert); // Import X values for solve
00381     X1 = (Epetra_Vector *) OverlapX_;
00382     Y1 = (Epetra_Vector *) OverlapY_; // Set pointers for X1 and Y1 to point to overlap space
00383   }
00384 
00385   Epetra_Flops * counter = this->GetFlopCounter();
00386   if (counter!=0) {
00387     L_->SetFlopCounter(*counter);
00388     Y1->SetFlopCounter(*counter);
00389     U_->SetFlopCounter(*counter);
00390   }
00391 
00392   if (!Trans) {
00393 
00394     L_->Solve(Lower, Trans, UnitDiagonal, *X1, *Y1);
00395     Y1->Multiply(1.0, *D_, *Y1, 0.0); // y = D*y (D_ has inverse of diagonal)
00396     U_->Solve(Upper, Trans, UnitDiagonal, *Y1, *Y1); // Solve Uy = y
00397   }
00398   else
00399     {
00400       U_->Solve(Upper, Trans, UnitDiagonal, *X1, *Y1); // Solve Uy = y
00401       Y1->Multiply(1.0, *D_, *Y1, 0.0); // y = D*y (D_ has inverse of diagonal)
00402       L_->Solve(Lower, Trans, UnitDiagonal, *Y1, *Y1);
00403       
00404     } 
00405 
00406   // Export computed Y values as directed
00407   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal())
00408     Y.Export(*OverlapY_,*Graph_.OverlapImporter(), OverlapMode_);
00409   return(0);
00410 }
00411 
00412 
00413 //=============================================================================
00414 int Ifpack_CrsRick::Solve(bool Trans, const Epetra_MultiVector& X, 
00415         Epetra_MultiVector& Y) const {
00416 //
00417 // This function finds Y such that LDU Y = X or U(trans) D L(trans) Y = X for multiple RHS
00418 //
00419 
00420   if (X.NumVectors()!=Y.NumVectors()) EPETRA_CHK_ERR(-1); // Return error: X and Y not the same size
00421 
00422   bool Upper = true;
00423   bool Lower = false;
00424   bool UnitDiagonal = true;
00425 
00426   Epetra_MultiVector * X1 = (Epetra_MultiVector *) &X;
00427   Epetra_MultiVector * Y1 = (Epetra_MultiVector *) &Y;
00428 
00429   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal()) {
00430     // Make sure the number of vectors in the multivector is the same as before.
00431     if (OverlapX_!=0) {
00432       if (OverlapX_->NumVectors()!=X.NumVectors()) {
00433   delete OverlapX_; OverlapX_ = 0;
00434   delete OverlapY_; OverlapY_ = 0;
00435       }
00436     }
00437     if (OverlapX_==0) { // Need to allocate space for overlap X and Y
00438       OverlapX_ = new Epetra_MultiVector(Graph_.OverlapGraph()->RowMap(), X.NumVectors());
00439       OverlapY_ = new Epetra_MultiVector(Graph_.OverlapGraph()->RowMap(), Y.NumVectors());
00440     }
00441     OverlapX_->Import(X,*Graph_.OverlapImporter(), Insert); // Import X values for solve
00442     X1 = OverlapX_;
00443     Y1 = OverlapY_; // Set pointers for X1 and Y1 to point to overlap space
00444   }
00445 
00446   Epetra_Flops * counter = this->GetFlopCounter();
00447   if (counter!=0) {
00448     L_->SetFlopCounter(*counter);
00449     Y1->SetFlopCounter(*counter);
00450     U_->SetFlopCounter(*counter);
00451   }
00452 
00453   if (!Trans) {
00454 
00455     L_->Solve(Lower, Trans, UnitDiagonal, *X1, *Y1);
00456     Y1->Multiply(1.0, *D_, *Y1, 0.0); // y = D*y (D_ has inverse of diagonal)
00457     U_->Solve(Upper, Trans, UnitDiagonal, *Y1, *Y1); // Solve Uy = y
00458   }
00459   else
00460     {
00461       U_->Solve(Upper, Trans, UnitDiagonal, *X1, *Y1); // Solve Uy = y
00462       Y1->Multiply(1.0, *D_, *Y1, 0.0); // y = D*y (D_ has inverse of diagonal)
00463       L_->Solve(Lower, Trans, UnitDiagonal, *Y1, *Y1);
00464       
00465     } 
00466 
00467   // Export computed Y values as directed
00468   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal())
00469     Y.Export(*OverlapY_,*Graph_.OverlapImporter(), OverlapMode_);
00470   return(0);
00471 }
00472 //=============================================================================
00473 int Ifpack_CrsRick::Multiply(bool Trans, const Epetra_MultiVector& X, 
00474         Epetra_MultiVector& Y) const {
00475 //
00476 // This function finds X such that LDU Y = X or U(trans) D L(trans) Y = X for multiple RHS
00477 //
00478 
00479   if (X.NumVectors()!=Y.NumVectors()) EPETRA_CHK_ERR(-1); // Return error: X and Y not the same size
00480 
00481   bool Upper = true;
00482   bool Lower = false;
00483   bool UnitDiagonal = true;
00484 
00485   Epetra_MultiVector * X1 = (Epetra_MultiVector *) &X;
00486   Epetra_MultiVector * Y1 = (Epetra_MultiVector *) &Y;
00487 
00488   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal()) {
00489     // Make sure the number of vectors in the multivector is the same as before.
00490     if (OverlapX_!=0) {
00491       if (OverlapX_->NumVectors()!=X.NumVectors()) {
00492   delete OverlapX_; OverlapX_ = 0;
00493   delete OverlapY_; OverlapY_ = 0;
00494       }
00495     }
00496     if (OverlapX_==0) { // Need to allocate space for overlap X and Y
00497       OverlapX_ = new Epetra_MultiVector(Graph_.OverlapGraph()->RowMap(), X.NumVectors());
00498       OverlapY_ = new Epetra_MultiVector(Graph_.OverlapGraph()->RowMap(), Y.NumVectors());
00499     }
00500     OverlapX_->Import(X,*Graph_.OverlapImporter(), Insert); // Import X values for solve
00501     X1 = OverlapX_;
00502     Y1 = OverlapY_; // Set pointers for X1 and Y1 to point to overlap space
00503   }
00504 
00505   Epetra_Flops * counter = this->GetFlopCounter();
00506   if (counter!=0) {
00507     L_->SetFlopCounter(*counter);
00508     Y1->SetFlopCounter(*counter);
00509     U_->SetFlopCounter(*counter);
00510   }
00511 
00512   if (Trans) {
00513 
00514     L_->Multiply(Trans, *X1, *Y1);
00515     Y1->Update(1.0, *X1, 1.0); // Y1 = Y1 + X1 (account for implicit unit diagonal)
00516     Y1->ReciprocalMultiply(1.0, *D_, *Y1, 0.0); // y = D*y (D_ has inverse of diagonal)
00517     Epetra_MultiVector Y1temp(*Y1); // Need a temp copy of Y1
00518     U_->Multiply(Trans, Y1temp, *Y1);
00519     Y1->Update(1.0, Y1temp, 1.0); // (account for implicit unit diagonal)
00520   }
00521   else {
00522     U_->Multiply(Trans, *X1, *Y1); // 
00523     Y1->Update(1.0, *X1, 1.0); // Y1 = Y1 + X1 (account for implicit unit diagonal)
00524     Y1->ReciprocalMultiply(1.0, *D_, *Y1, 0.0); // y = D*y (D_ has inverse of diagonal)
00525     Epetra_MultiVector Y1temp(*Y1); // Need a temp copy of Y1
00526     L_->Multiply(Trans, Y1temp, *Y1);
00527     Y1->Update(1.0, Y1temp, 1.0); // (account for implicit unit diagonal)
00528     } 
00529 
00530   // Export computed Y values as directed
00531   if (Graph_.LevelOverlap()>0 && Graph_.L_Graph().DomainMap().DistributedGlobal())
00532     Y.Export(*OverlapY_,*Graph_.OverlapImporter(), OverlapMode_);
00533   return(0);
00534 }
00535 //=============================================================================
00536 int Ifpack_CrsRick::Condest(bool Trans, double & ConditionNumberEstimate) const {
00537 
00538   if (Condest_>=0.0) {
00539     ConditionNumberEstimate = Condest_;
00540     return(0);
00541   }
00542   // Create a vector with all values equal to one
00543   Epetra_Vector Ones(A_.RowMap());
00544   Epetra_Vector OnesResult(Ones);
00545   Ones.PutScalar(1.0);
00546 
00547   EPETRA_CHK_ERR(Solve(Trans, Ones, OnesResult)); // Compute the effect of the solve on the vector of ones
00548   EPETRA_CHK_ERR(OnesResult.Abs(OnesResult)); // Make all values non-negative
00549   EPETRA_CHK_ERR(OnesResult.MaxValue(&ConditionNumberEstimate)); // Get the maximum value across all processors
00550   Condest_ = ConditionNumberEstimate; // Save value for possible later calls
00551   return(0);
00552 }
00553 //=============================================================================
00554 // Non-member functions
00555 
00556 ostream& operator << (ostream& os, const Ifpack_CrsRick& A)
00557 {
00558 /*  Epetra_fmtflags olda = os.setf(ios::right,ios::adjustfield);
00559   Epetra_fmtflags oldf = os.setf(ios::scientific,ios::floatfield);
00560   int oldp = os.precision(12); */
00561   int LevelFill = A.Graph().LevelFill();
00562   int LevelOverlap = A.Graph().LevelOverlap();
00563   Epetra_CrsMatrix & L = (Epetra_CrsMatrix &) A.L();
00564   Epetra_CrsMatrix & U = (Epetra_CrsMatrix &) A.U();
00565   Epetra_Vector & D = (Epetra_Vector &) A.D();
00566 
00567   os.width(14);
00568   os << endl;
00569   os <<  "     Level of Fill = "; os << LevelFill;
00570   os << endl;
00571   os.width(14);
00572   os <<  "     Level of Overlap = "; os << LevelOverlap;
00573   os << endl;
00574 
00575   os.width(14);
00576   os <<  "     Lower Triangle = ";
00577   os << endl;
00578   os << L; // Let Epetra_CrsMatrix handle the rest.
00579   os << endl;
00580 
00581   os.width(14);
00582   os <<  "     Inverse of Diagonal = ";
00583   os << endl;
00584   os << D; // Let Epetra_Vector handle the rest.
00585   os << endl;
00586 
00587   os.width(14);
00588   os <<  "     Upper Triangle = ";
00589   os << endl;
00590   os << U; // Let Epetra_CrsMatrix handle the rest.
00591   os << endl;
00592  
00593   // Reset os flags
00594 
00595 /*  os.setf(olda,ios::adjustfield);
00596   os.setf(oldf,ios::floatfield);
00597   os.precision(oldp); */
00598 
00599   return os;
00600 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines