Epetra Package Browser (Single Doxygen Collection) Development
example/CrsSingletonFilterExample/cxx_main.cpp
Go to the documentation of this file.
00001 //@HEADER
00002 // ************************************************************************
00003 //
00004 //               Epetra: Linear Algebra Services Package
00005 //                 Copyright 2011 Sandia Corporation
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ************************************************************************
00040 //@HEADER
00041 
00042 #include "Epetra_Comm.h"
00043 #include "Epetra_Map.h"
00044 #include "Epetra_Time.h"
00045 #include "Epetra_BlockMap.h"
00046 #include "Epetra_MultiVector.h"
00047 #include "Epetra_Vector.h"
00048 #include "Epetra_Export.h"
00049 #include "Epetra_LinearProblem.h"
00050 #include "Epetra_CrsSingletonFilter.h"
00051 
00052 #include "Epetra_VbrMatrix.h"
00053 #include "Epetra_CrsMatrix.h"
00054 #ifdef EPETRA_MPI
00055 #include "mpi.h"
00056 #include "Epetra_MpiComm.h"
00057 #else
00058 #include "Epetra_SerialComm.h"
00059 #endif
00060 #include "Trilinos_Util.h"
00061 #include "Ifpack_IlukGraph.h"
00062 #include "Ifpack_CrsRiluk.h"
00063 
00064 void BiCGSTAB(Epetra_CrsMatrix &A, Epetra_Vector &x, Epetra_Vector &b,
00065         Ifpack_CrsRiluk *M,
00066         int Maxiter, double Tolerance,
00067         double *residual, bool verbose);
00068 int Statistics(const Epetra_CrsSingletonFilter & filter);
00069 int main(int argc, char *argv[]) {
00070 
00071 #ifdef EPETRA_MPI
00072   MPI_Init(&argc,&argv);
00073   Epetra_MpiComm Comm (MPI_COMM_WORLD);
00074 #else
00075   Epetra_SerialComm Comm;
00076 #endif
00077 
00078   cout << Comm << endl;
00079 
00080   int MyPID = Comm.MyPID();
00081 
00082   bool verbose = false;
00083   bool verbose1 = true;
00084   if (MyPID==0) verbose = true;
00085 
00086   if(argc < 2 && verbose) {
00087     cerr << "Usage: " << argv[0]
00088    << " HB_filename [level_fill [level_overlap [absolute_threshold [ relative_threshold]]]]" << endl
00089    << "where:" << endl
00090    << "HB_filename        - filename and path of a Harwell-Boeing data set" << endl
00091    << "level_fill         - The amount of fill to use for ILU(k) preconditioner (default 0)" << endl
00092    << "level_overlap      - The amount of overlap used for overlapping Schwarz subdomains (default 0)" << endl
00093    << "absolute_threshold - The minimum value to place on the diagonal prior to factorization (default 0.0)" << endl
00094    << "relative_threshold - The relative amount to perturb the diagonal prior to factorization (default 1.0)" << endl << endl
00095    << "To specify a non-default value for one of these parameters, you must specify all" << endl
00096    << " preceding values but not any subsequent parameters. Example:" << endl
00097    << "ifpackHpcSerialMsr.exe mymatrix.hpc 1  - loads mymatrix.hpc, uses level fill of one, all other values are defaults" << endl
00098    << endl;
00099     return(1);
00100 
00101   }
00102 
00103   // Uncomment the next three lines to debug in mpi mode
00104   //int tmp;
00105   //if (MyPID==0) cin >> tmp;
00106   //Comm.Barrier();
00107 
00108   Epetra_Map * readMap;
00109   Epetra_CrsMatrix * readA;
00110   Epetra_Vector * readx;
00111   Epetra_Vector * readb;
00112   Epetra_Vector * readxexact;
00113 
00114   // Call routine to read in HB problem
00115   Trilinos_Util_ReadHb2Epetra(argv[1], Comm, readMap, readA, readx, readb, readxexact);
00116 
00117   // Create uniform distributed map
00118   Epetra_Map map(readMap->NumGlobalElements(), 0, Comm);
00119 
00120   // Create Exporter to distribute read-in matrix and vectors
00121 
00122   Epetra_Export exporter(*readMap, map);
00123   Epetra_CrsMatrix A(Copy, map, 0);
00124   Epetra_Vector x(map);
00125   Epetra_Vector b(map);
00126   Epetra_Vector xexact(map);
00127 
00128   Epetra_Time FillTimer(Comm);
00129   x.Export(*readx, exporter, Add);
00130   b.Export(*readb, exporter, Add);
00131   xexact.Export(*readxexact, exporter, Add);
00132   Comm.Barrier();
00133   double vectorRedistributeTime = FillTimer.ElapsedTime();
00134   A.Export(*readA, exporter, Add);
00135   Comm.Barrier();
00136   double matrixRedistributeTime = FillTimer.ElapsedTime() - vectorRedistributeTime;
00137   assert(A.FillComplete()==0);
00138   Comm.Barrier();
00139   double fillCompleteTime = FillTimer.ElapsedTime() - matrixRedistributeTime;
00140   if (Comm.MyPID()==0)  {
00141     cout << "\n\n****************************************************" << endl;
00142     cout << "\n Vector redistribute  time (sec) = " << vectorRedistributeTime<< endl;
00143     cout << "    Matrix redistribute time (sec) = " << matrixRedistributeTime << endl;
00144     cout << "    Transform to Local  time (sec) = " << fillCompleteTime << endl<< endl;
00145   }
00146   Epetra_Vector tmp1(*readMap);
00147   Epetra_Vector tmp2(map);
00148   readA->Multiply(false, *readxexact, tmp1);
00149 
00150   A.Multiply(false, xexact, tmp2);
00151   double residual;
00152   tmp1.Norm2(&residual);
00153   if (verbose) cout << "Norm of Ax from file            = " << residual << endl;
00154   tmp2.Norm2(&residual);
00155   if (verbose) cout << "Norm of Ax after redistribution = " << residual << endl << endl << endl;
00156 
00157   //cout << "A from file = " << *readA << endl << endl << endl;
00158 
00159   //cout << "A after dist = " << A << endl << endl << endl;
00160 
00161   delete readA;
00162   delete readx;
00163   delete readb;
00164   delete readxexact;
00165   delete readMap;
00166 
00167   Comm.Barrier();
00168 
00169   bool smallProblem = false;
00170   if (A.RowMap().NumGlobalElements()<100) smallProblem = true;
00171 
00172   if (smallProblem)
00173     cout << "Original Matrix = " << endl << A   << endl;
00174 
00175   x.PutScalar(0.0);
00176 
00177   Epetra_LinearProblem FullProblem(&A, &x, &b);
00178   double normb, norma;
00179   b.NormInf(&normb);
00180   norma = A.NormInf();
00181   if (verbose)
00182     cout << "Inf norm of Original Matrix = " << norma << endl
00183    << "Inf norm of Original RHS    = " << normb << endl;
00184 
00185   Epetra_Time ReductionTimer(Comm);
00186   Epetra_CrsSingletonFilter SingletonFilter;
00187   Comm.Barrier();
00188   double reduceInitTime = ReductionTimer.ElapsedTime();
00189   SingletonFilter.Analyze(&A);
00190   Comm.Barrier();
00191   double reduceAnalyzeTime = ReductionTimer.ElapsedTime() - reduceInitTime;
00192 
00193   if (SingletonFilter.SingletonsDetected())
00194     cout << "Singletons found" << endl;
00195   else {
00196     cout << "Singletons not found" << endl;
00197     exit(1);
00198   }
00199   SingletonFilter.ConstructReducedProblem(&FullProblem);
00200   Comm.Barrier();
00201   double reduceConstructTime = ReductionTimer.ElapsedTime() - reduceInitTime;
00202 
00203   double totalReduceTime = ReductionTimer.ElapsedTime();
00204 
00205   if (verbose)
00206     cout << "\n\n****************************************************" << endl
00207    << "    Reduction init  time (sec)           = " << reduceInitTime<< endl
00208    << "    Reduction Analyze time (sec)         = " << reduceAnalyzeTime << endl
00209    << "    Construct Reduced Problem time (sec) = " << reduceConstructTime << endl
00210    << "    Reduction Total time (sec)           = " << totalReduceTime << endl<< endl;
00211 
00212   Statistics(SingletonFilter);
00213 
00214   Epetra_LinearProblem * ReducedProblem = SingletonFilter.ReducedProblem();
00215 
00216   Epetra_CrsMatrix * Ap = dynamic_cast<Epetra_CrsMatrix *>(ReducedProblem->GetMatrix());
00217   Epetra_Vector * bp = (*ReducedProblem->GetRHS())(0);
00218   Epetra_Vector * xp = (*ReducedProblem->GetLHS())(0);
00219 
00220 
00221   if (smallProblem)
00222     cout << " Reduced Matrix = " << endl << *Ap << endl
00223    << " LHS before sol = " << endl << *xp << endl
00224    << " RHS            = " << endl << *bp << endl;
00225 
00226   // Construct ILU preconditioner
00227 
00228   double elapsed_time, total_flops, MFLOPs;
00229   Epetra_Time timer(Comm);
00230 
00231   int LevelFill = 0;
00232   if (argc > 2)  LevelFill = atoi(argv[2]);
00233   if (verbose) cout << "Using Level Fill = " << LevelFill << endl;
00234   int Overlap = 0;
00235   if (argc > 3) Overlap = atoi(argv[3]);
00236   if (verbose) cout << "Using Level Overlap = " << Overlap << endl;
00237   double Athresh = 0.0;
00238   if (argc > 4) Athresh = atof(argv[4]);
00239   if (verbose) cout << "Using Absolute Threshold Value of = " << Athresh << endl;
00240 
00241   double Rthresh = 1.0;
00242   if (argc > 5) Rthresh = atof(argv[5]);
00243   if (verbose) cout << "Using Relative Threshold Value of = " << Rthresh << endl;
00244 
00245   Ifpack_IlukGraph * IlukGraph = 0;
00246   Ifpack_CrsRiluk * ILUK = 0;
00247 
00248   if (LevelFill>-1) {
00249     elapsed_time = timer.ElapsedTime();
00250     IlukGraph = new Ifpack_IlukGraph(Ap->Graph(), LevelFill, Overlap);
00251     assert(IlukGraph->ConstructFilledGraph()==0);
00252     elapsed_time = timer.ElapsedTime() - elapsed_time;
00253     if (verbose) cout << "Time to construct ILUK graph = " << elapsed_time << endl;
00254 
00255 
00256     Epetra_Flops fact_counter;
00257 
00258     elapsed_time = timer.ElapsedTime();
00259     ILUK = new Ifpack_CrsRiluk(*IlukGraph);
00260     ILUK->SetFlopCounter(fact_counter);
00261     ILUK->SetAbsoluteThreshold(Athresh);
00262     ILUK->SetRelativeThreshold(Rthresh);
00263     //assert(ILUK->InitValues()==0);
00264     int initerr = ILUK->InitValues(*Ap);
00265     if (initerr!=0) {
00266       cout << endl << Comm << endl << "  InitValues error = " << initerr;
00267       if (initerr==1) cout << "  Zero diagonal found, warning error only";
00268       cout << endl << endl;
00269     }
00270     assert(ILUK->Factor()==0);
00271     elapsed_time = timer.ElapsedTime() - elapsed_time;
00272     total_flops = ILUK->Flops();
00273     MFLOPs = total_flops/elapsed_time/1000000.0;
00274     if (verbose) cout << "Time to compute preconditioner values = "
00275         << elapsed_time << endl
00276         << "MFLOPS for Factorization = " << MFLOPs << endl;
00277     //cout << *ILUK << endl;
00278   double Condest;
00279   ILUK->Condest(false, Condest);
00280 
00281   if (verbose) cout << "Condition number estimate for this preconditioner = " << Condest << endl;
00282   }
00283   int Maxiter = 100;
00284   double Tolerance = 1.0E-8;
00285 
00286   Epetra_Flops counter;
00287   Ap->SetFlopCounter(counter);
00288   xp->SetFlopCounter(*Ap);
00289   bp->SetFlopCounter(*Ap);
00290   if (ILUK!=0) ILUK->SetFlopCounter(*Ap);
00291 
00292   elapsed_time = timer.ElapsedTime();
00293 
00294   double normreducedb, normreduceda;
00295   bp->NormInf(&normreducedb);
00296   normreduceda = Ap->NormInf();
00297   if (verbose)
00298     cout << "Inf norm of Reduced Matrix = " << normreduceda << endl
00299    << "Inf norm of Reduced RHS    = " << normreducedb << endl;
00300 
00301   BiCGSTAB(*Ap, *xp, *bp, ILUK, Maxiter, Tolerance, &residual, verbose);
00302 
00303   elapsed_time = timer.ElapsedTime() - elapsed_time;
00304   total_flops = counter.Flops();
00305   MFLOPs = total_flops/elapsed_time/1000000.0;
00306   if (verbose) cout << "Time to compute solution = "
00307         << elapsed_time << endl
00308         << "Number of operations in solve = " << total_flops << endl
00309         << "MFLOPS for Solve = " << MFLOPs<< endl << endl;
00310 
00311   SingletonFilter.ComputeFullSolution();
00312 
00313   if (smallProblem)
00314   cout << " Reduced LHS after sol = " << endl << *xp << endl
00315        << " Full    LHS after sol = " << endl << x << endl
00316        << " Full  Exact LHS         = " << endl << xexact << endl;
00317 
00318   Epetra_Vector resid(x);
00319 
00320   resid.Update(1.0, x, -1.0, xexact, 0.0); // resid = xcomp - xexact
00321 
00322   resid.Norm2(&residual);
00323   double normx, normxexact;
00324   x.Norm2(&normx);
00325   xexact.Norm2(&normxexact);
00326 
00327   if (verbose)
00328     cout << "2-norm of computed solution                               = " << normx << endl
00329    << "2-norm of exact solution                                  = " << normxexact << endl
00330    << "2-norm of difference between computed and exact solution  = " << residual << endl;
00331 
00332   if (verbose1 && residual>1.0e-5) {
00333     if (verbose)
00334       cout << "Difference between computed and exact solution appears large..." << endl
00335      << "Computing norm of A times this difference.  If this norm is small, then matrix is singular"
00336      << endl;
00337     Epetra_Vector bdiff(b);
00338     assert(A.Multiply(false, resid, bdiff)==0);
00339     assert(bdiff.Norm2(&residual)==0);
00340     if (verbose)
00341       cout << "2-norm of A times difference between computed and exact solution  = " << residual << endl;
00342 
00343   }
00344   if (verbose)
00345     cout << "********************************************************" << endl
00346    << "              Solving again with 2*Ax=2*b" << endl
00347    << "********************************************************" << endl;
00348 
00349   A.Scale(1.0); // A = 2*A
00350   b.Scale(1.0); // b = 2*b
00351   x.PutScalar(0.0);
00352   b.NormInf(&normb);
00353   norma = A.NormInf();
00354   if (verbose)
00355     cout << "Inf norm of Original Matrix = " << norma << endl
00356    << "Inf norm of Original RHS    = " << normb << endl;
00357   double updateReducedProblemTime = ReductionTimer.ElapsedTime();
00358   SingletonFilter.UpdateReducedProblem(&FullProblem);
00359   Comm.Barrier();
00360   updateReducedProblemTime = ReductionTimer.ElapsedTime() - updateReducedProblemTime;
00361   if (verbose)
00362     cout << "\n\n****************************************************" << endl
00363    << "    Update Reduced Problem time (sec)           = " << updateReducedProblemTime<< endl
00364    << "****************************************************" << endl;
00365   Statistics(SingletonFilter);
00366 
00367   if (LevelFill>-1) {
00368 
00369     Epetra_Flops fact_counter;
00370 
00371     elapsed_time = timer.ElapsedTime();
00372 
00373     int initerr = ILUK->InitValues(*Ap);
00374     if (initerr!=0) {
00375       cout << endl << Comm << endl << "  InitValues error = " << initerr;
00376       if (initerr==1) cout << "  Zero diagonal found, warning error only";
00377       cout << endl << endl;
00378     }
00379     assert(ILUK->Factor()==0);
00380     elapsed_time = timer.ElapsedTime() - elapsed_time;
00381     total_flops = ILUK->Flops();
00382     MFLOPs = total_flops/elapsed_time/1000000.0;
00383     if (verbose) cout << "Time to compute preconditioner values = "
00384         << elapsed_time << endl
00385         << "MFLOPS for Factorization = " << MFLOPs << endl;
00386     double Condest;
00387     ILUK->Condest(false, Condest);
00388 
00389     if (verbose) cout << "Condition number estimate for this preconditioner = " << Condest << endl;
00390   }
00391   bp->NormInf(&normreducedb);
00392   normreduceda = Ap->NormInf();
00393   if (verbose)
00394     cout << "Inf norm of Reduced Matrix = " << normreduceda << endl
00395    << "Inf norm of Reduced RHS    = " << normreducedb << endl;
00396 
00397   BiCGSTAB(*Ap, *xp, *bp, ILUK, Maxiter, Tolerance, &residual, verbose);
00398 
00399   elapsed_time = timer.ElapsedTime() - elapsed_time;
00400   total_flops = counter.Flops();
00401   MFLOPs = total_flops/elapsed_time/1000000.0;
00402   if (verbose) cout << "Time to compute solution = "
00403         << elapsed_time << endl
00404         << "Number of operations in solve = " << total_flops << endl
00405         << "MFLOPS for Solve = " << MFLOPs<< endl << endl;
00406 
00407   SingletonFilter.ComputeFullSolution();
00408 
00409   if (smallProblem)
00410   cout << " Reduced LHS after sol = " << endl << *xp << endl
00411        << " Full    LHS after sol = " << endl << x << endl
00412        << " Full  Exact LHS         = " << endl << xexact << endl;
00413 
00414   resid.Update(1.0, x, -1.0, xexact, 0.0); // resid = xcomp - xexact
00415 
00416   resid.Norm2(&residual);
00417   x.Norm2(&normx);
00418   xexact.Norm2(&normxexact);
00419 
00420   if (verbose)
00421     cout << "2-norm of computed solution                               = " << normx << endl
00422    << "2-norm of exact solution                                  = " << normxexact << endl
00423    << "2-norm of difference between computed and exact solution  = " << residual << endl;
00424 
00425   if (verbose1 && residual>1.0e-5) {
00426     if (verbose)
00427       cout << "Difference between computed and exact solution appears large..." << endl
00428      << "Computing norm of A times this difference.  If this norm is small, then matrix is singular"
00429      << endl;
00430     Epetra_Vector bdiff(b);
00431     assert(A.Multiply(false, resid, bdiff)==0);
00432     assert(bdiff.Norm2(&residual)==0);
00433     if (verbose)
00434       cout << "2-norm of A times difference between computed and exact solution  = " << residual << endl;
00435 
00436   }
00437 
00438 
00439 
00440   if (ILUK!=0) delete ILUK;
00441   if (IlukGraph!=0) delete IlukGraph;
00442         
00443 #ifdef EPETRA_MPI
00444   MPI_Finalize() ;
00445 #endif
00446 
00447 return 0 ;
00448 }
00449 void BiCGSTAB(Epetra_CrsMatrix &A,
00450         Epetra_Vector &x,
00451         Epetra_Vector &b,
00452         Ifpack_CrsRiluk *M,
00453         int Maxiter,
00454         double Tolerance,
00455         double *residual, bool verbose) {
00456 
00457   // Allocate vectors needed for iterations
00458   Epetra_Vector phat(x.Map()); phat.SetFlopCounter(x);
00459   Epetra_Vector p(x.Map()); p.SetFlopCounter(x);
00460   Epetra_Vector shat(x.Map()); shat.SetFlopCounter(x);
00461   Epetra_Vector s(x.Map()); s.SetFlopCounter(x);
00462   Epetra_Vector r(b.Map()); r.SetFlopCounter(x);
00463   Epetra_Vector rtilde(x.Map()); rtilde.PutScalar(1.0); rtilde.SetFlopCounter(x);
00464   Epetra_Vector v(x.Map()); v.SetFlopCounter(x);
00465 
00466 
00467   A.Multiply(false, x, r); // r = A*x
00468 
00469   r.Update(1.0, b, -1.0); // r = b - A*x
00470 
00471   double r_norm, b_norm, scaled_r_norm, rhon, rhonm1 = 1.0;
00472   double alpha = 1.0, omega = 1.0, sigma;
00473   double omega_num, omega_den;
00474   r.Norm2(&r_norm);
00475   b.Norm2(&b_norm);
00476   scaled_r_norm = r_norm/b_norm;
00477   r.Dot(rtilde,&rhon);
00478 
00479   if (verbose) cout << "Initial residual = " << r_norm
00480         << " Scaled residual = " << scaled_r_norm << endl;
00481 
00482 
00483   for (int i=0; i<Maxiter; i++) { // Main iteration loop
00484 
00485     double beta = (rhon/rhonm1) * (alpha/omega);
00486     rhonm1 = rhon;
00487 
00488     /* p    = r + beta*(p - omega*v)       */
00489     /* phat = M^-1 p                       */
00490     /* v    = A phat                       */
00491 
00492     double dtemp = - beta*omega;
00493 
00494     p.Update(1.0, r, dtemp, v, beta);
00495     if (M==0)
00496       phat.Scale(1.0, p);
00497     else
00498       M->Solve(false, p, phat);
00499     A.Multiply(false, phat, v);
00500 
00501 
00502     rtilde.Dot(v,&sigma);
00503     alpha = rhon/sigma;
00504 
00505     /* s = r - alpha*v                     */
00506     /* shat = M^-1 s                       */
00507     /* r = A shat (r is a tmp here for t ) */
00508 
00509     s.Update(-alpha, v, 1.0, r, 0.0);
00510     if (M==0)
00511       shat.Scale(1.0, s);
00512     else
00513       M->Solve(false, s, shat);
00514     A.Multiply(false, shat, r);
00515 
00516     r.Dot(s, &omega_num);
00517     r.Dot(r, &omega_den);
00518     omega = omega_num / omega_den;
00519 
00520     /* x = x + alpha*phat + omega*shat */
00521     /* r = s - omega*r */
00522 
00523     x.Update(alpha, phat, omega, shat, 1.0);
00524     r.Update(1.0, s, -omega);
00525 
00526     r.Norm2(&r_norm);
00527     scaled_r_norm = r_norm/b_norm;
00528     r.Dot(rtilde,&rhon);
00529 
00530     if (verbose) cout << "Iter "<< i << " residual = " << r_norm
00531           << " Scaled residual = " << scaled_r_norm << endl;
00532 
00533     if (r_norm < Tolerance) break;
00534   }
00535   return;
00536 }
00537 //==============================================================================
00538 int Statistics(const Epetra_CrsSingletonFilter & filter) {
00539 
00540   // Create local variables of some of the stats we need because we don't know if the
00541   // method calls are collective or not for the particular implementation of the Row Matrix
00542   int fmaxentries;
00543   int maxentries = filter.FullMatrix()->MaxNumEntries();
00544   filter.FullMatrix()->RowMatrixRowMap().Comm().SumAll(&maxentries, &fmaxentries, 1);
00545   int fnrows = filter.FullMatrix()->NumGlobalRows();
00546   int fnnz = filter.FullMatrix()->NumGlobalNonzeros();
00547   if (filter.FullMatrix()->RowMatrixRowMap().Comm().MyPID()!=0) return(0);
00548 
00549   cout << "Full System characteristics:" << endl << endl
00550        << "  Dimension                             = " << fnrows << endl
00551        << "  Number of nonzeros                    = " <<fnnz << endl
00552        << "  Maximum Number of Row Entries         = " << fmaxentries << endl << endl
00553        << "Reduced System characteristics:" << endl << endl
00554        << "  Dimension                             = " << filter.ReducedMatrix()->NumGlobalRows() << endl
00555        << "  Number of nonzeros                    = " << filter.ReducedMatrix()->NumGlobalNonzeros() << endl
00556        << "  Maximum Number of Row Entries         = " << filter.ReducedMatrix()->GlobalMaxNumEntries() << endl << endl
00557        << "Singleton information: " << endl
00558        << "  Total number of singletons            = " << filter.NumSingletons() << endl
00559        << "  Number of rows with single entries    = " << filter.NumRowSingletons() << endl
00560        << "  Number of columns with single entries " << endl
00561        << "    (that were not already counted as " << endl
00562        << "     row singletons)                    = " << filter.NumColSingletons() << endl << endl
00563        << "Ratios: " << endl
00564        << "  Percent reduction in dimension        = " << filter.RatioOfDimensions()*100.0 << endl
00565        << "  Percent reduction in nonzero count    = " << filter.RatioOfNonzeros()*100.0 << endl << endl;
00566 
00567   return(0);
00568 
00569 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines