Epetra Package Browser (Single Doxygen Collection) Development
example/InverseIteration/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 <stdio.h>
00043 #include <stdlib.h>
00044 #include <ctype.h>
00045 #include <assert.h>
00046 #include <string.h>
00047 #include <math.h>
00048 #include "Epetra_Comm.h"
00049 #include "Epetra_Map.h"
00050 #include "Epetra_Time.h"
00051 #include "Epetra_BlockMap.h"
00052 #include "Epetra_MultiVector.h"
00053 #include "Epetra_Vector.h"
00054 #include "Epetra_CrsMatrix.h"
00055 #ifdef EPETRA_MPI
00056 #include "mpi.h"
00057 #include "Epetra_MpiComm.h"
00058 #else
00059 #include "Epetra_SerialComm.h"
00060 #endif
00061 #include "Trilinos_Util.h"
00062 #include "Ifpack_IlukGraph.h"
00063 #include "Ifpack_CrsRiluk.h"
00064 
00065 // *************************************************************
00066 // Function prototypes
00067 // *************************************************************
00068 
00069 int invIteration(Epetra_CrsMatrix& A, double &lambda, bool verbose);
00070 int applyInverseSetup(Epetra_CrsMatrix &A, Ifpack_CrsRiluk * & M);
00071 int applyInverse(Epetra_CrsMatrix &A, Epetra_Vector & x, Epetra_Vector & b, Ifpack_CrsRiluk * M,
00072      bool verbose);
00073 int applyInverseDestroy(Ifpack_CrsRiluk * M);
00074 void BiCGSTAB(Epetra_CrsMatrix &A,
00075         Epetra_Vector &x,
00076         Epetra_Vector &b,
00077         Ifpack_CrsRiluk *M,
00078         int Maxiter,
00079         double Tolerance, bool verbose);
00080 
00081 // *************************************************************
00082 // main program - This benchmark code reads a Harwell-Boeing data
00083 //                set and finds the minimal eigenvalue of the matrix
00084 //                using inverse iteration.
00085 // *************************************************************
00086 int main(int argc, char *argv[]) {
00087 
00088 #ifdef EPETRA_MPI
00089   MPI_Init(&argc,&argv);
00090   Epetra_MpiComm Comm (MPI_COMM_WORLD);
00091 #else
00092   Epetra_SerialComm Comm;
00093 #endif
00094 
00095   cout << Comm << endl;
00096 
00097   int MyPID = Comm.MyPID();
00098 
00099   bool verbose = false;
00100   if (MyPID==0) verbose = true; // Print out detailed results (turn off for best performance)
00101 
00102   if(argc != 2) {
00103     if (verbose) cerr << "Usage: " << argv[0] << " HB_data_file" << endl;
00104     exit(1); // Error
00105   }
00106 
00107   // Define pointers that will be set by HB read function
00108 
00109   Epetra_Map * readMap;
00110   Epetra_CrsMatrix * readA;
00111   Epetra_Vector * readx;
00112   Epetra_Vector * readb;
00113   Epetra_Vector * readxexact;
00114 
00115   // Call function to read in HB problem
00116   Trilinos_Util_ReadHb2Epetra(argv[1], Comm, readMap, readA, readx, readb, readxexact);
00117 
00118   // Not interested in x, b or xexact for an eigenvalue problem
00119   delete readx;
00120   delete readb;
00121   delete readxexact;
00122 
00123 #ifdef EPETRA_MPI // If running in parallel, we need to distribute matrix across all PEs.
00124 
00125   // Create uniform distributed map
00126   Epetra_Map map(readMap->NumGlobalElements(), 0, Comm);
00127 
00128   // Create Exporter to distribute read-in matrix and vectors
00129 
00130   Epetra_Export exporter(*readMap, map);
00131   Epetra_CrsMatrix A(Copy, map, 0);
00132 
00133   A.Export(*readA, exporter, Add);
00134   assert(A.FillComplete()==0);
00135 
00136   delete readA;
00137   delete readMap;
00138 
00139 #else // If not running in parallel, we do not need to distribute the matrix
00140   Epetra_CrsMatrix & A = *readA;
00141 #endif
00142 
00143   // Create flop counter to collect all FLOPS
00144   Epetra_Flops counter;
00145   A.SetFlopCounter(counter);
00146 
00147   double lambda = 0; // Minimal eigenvalue returned here
00148   // Call inverse iteration solver
00149   Epetra_Time timer(Comm);
00150   invIteration(A, lambda, verbose);
00151   double elapsedTime = timer.ElapsedTime();
00152   double totalFlops = counter.Flops();
00153   double MFLOPS = totalFlops/elapsedTime/1000000.0;
00154 
00155 
00156   cout << endl
00157        << "*************************************************" << endl
00158        << " Approximate smallest eigenvalue = " << lambda << endl
00159        << "    Total Time    = " << elapsedTime << endl
00160        << "    Total FLOPS   = " << totalFlops << endl
00161        << "    Total MFLOPS  = " << MFLOPS << endl
00162        << "*************************************************" << endl;
00163 
00164   // All done
00165 #ifdef EPETRA_MPI
00166   MPI_Finalize();
00167 #else
00168   delete readA;
00169   delete readMap;
00170 #endif
00171 
00172 return (0);
00173 }
00174 
00175 // The following functions are used to solve the problem:
00176 
00177 
00178 // *************************************************************
00179 // Computes the smallest eigenvalue of the given matrix A.
00180 // Returns result in lambda
00181 // *************************************************************
00182 
00183 int invIteration(Epetra_CrsMatrix& A, double &lambda, bool verbose) {
00184 
00185   Ifpack_CrsRiluk * M;
00186   applyInverseSetup(A, M);
00187   Epetra_Vector q(A.RowMap());
00188   Epetra_Vector z(A.RowMap());
00189   Epetra_Vector resid(A.RowMap());
00190 
00191   Epetra_Flops * counter = A.GetFlopCounter();
00192   if (counter!=0) {
00193     q.SetFlopCounter(A);
00194     z.SetFlopCounter(A);
00195     resid.SetFlopCounter(A);
00196   }
00197 
00198   // Fill z with random Numbers
00199   z.Random();
00200 
00201   // variable needed for iteration
00202   double normz, residual;
00203 
00204   int niters = 100;
00205   double tolerance = 1.0E-6;
00206   int ierr = 1;
00207 
00208   for (int iter = 0; iter < niters; iter++)
00209     {
00210       if (verbose)
00211   cout << endl
00212        << " ***** Performing step " << iter << " of inverse iteration ***** " << endl;
00213 
00214       z.Norm2(&normz); // Compute 2-norm of z
00215       q.Scale(1.0/normz, z);
00216       applyInverse(A, z, q, M, verbose); // Compute z such that Az = q
00217       q.Dot(z, &lambda); // Approximate maximum eigenvalue
00218       if (iter%10==0 || iter+1==niters)
00219   {
00220     resid.Update(1.0, z, -lambda, q, 0.0); // Compute A(inv)*q - lambda*q
00221     resid.Norm2(&residual);
00222     cout << endl
00223          << "***** Inverse Iteration Step " << iter+1 << endl
00224          << "  Lambda = " << 1.0/lambda << endl
00225          << "  Residual of A(inv)*q - lambda*q = "
00226          << residual << endl;
00227   }
00228       if (residual < tolerance) {
00229   ierr = 0;
00230   break;
00231       }
00232     }
00233   // lambda is the largest eigenvalue of A(inv).  1/lambda is smallest eigenvalue of A.
00234   lambda = 1.0/lambda;
00235 
00236   // Compute A*q - lambda*q explicitly
00237   A.Multiply(false, q, z);
00238   resid.Update(1.0, z, -lambda, q, 0.0); // Compute A*q - lambda*q
00239   resid.Norm2(&residual);
00240   cout << "  Explicitly computed residual of A*q - lambda*q = " << residual << endl;
00241   applyInverseDestroy(M);
00242   return(ierr);
00243 }
00244 
00245 // Performs any setup that can be done once to reduce the cost of the applyInverse function
00246 int applyInverseSetup(Epetra_CrsMatrix &A, Ifpack_CrsRiluk * & M) {
00247   int LevelFill = 4;
00248   int Overlap = 0;
00249   Ifpack_IlukGraph * IlukGraph = new Ifpack_IlukGraph(A.Graph(), LevelFill, Overlap);
00250     assert(IlukGraph->ConstructFilledGraph()==0);
00251     M = new Ifpack_CrsRiluk(A, *IlukGraph);
00252     M->SetFlopCounter(A);
00253     assert(M->InitValues()==0);
00254     assert(M->Factor()==0);
00255   return(0);
00256 }
00257 
00258 // *************************************************************
00259 // Solves a problem Ax = b, for a given A and b.
00260 // M is a preconditioner computed in applyInverseSetup.
00261 // *************************************************************
00262 
00263 int applyInverse(Epetra_CrsMatrix &A,
00264      Epetra_Vector & x, Epetra_Vector & b, Ifpack_CrsRiluk * M, bool verbose) {
00265 
00266   int Maxiter = 1000;
00267   double Tolerance = 1.0E-16;
00268   BiCGSTAB(A, x, b, M, Maxiter, Tolerance, verbose);
00269 
00270   return(0);
00271 }
00272 
00273 
00274 // *************************************************************
00275 // Releases any memory associated with the preconditioner M
00276 // *************************************************************
00277 
00278 int applyInverseDestroy(Ifpack_CrsRiluk * M) {
00279   Ifpack_IlukGraph & IlukGraph = (Ifpack_IlukGraph &) M->Graph();
00280   delete M;
00281   delete &IlukGraph;
00282   return(0);
00283 }
00284 
00285 // *************************************************************
00286 // Uses the Bi-CGSTAB iterative method to solve Ax = b
00287 // *************************************************************
00288 
00289 void BiCGSTAB(Epetra_CrsMatrix &A,
00290         Epetra_Vector &x,
00291         Epetra_Vector &b,
00292         Ifpack_CrsRiluk *M,
00293         int Maxiter,
00294         double Tolerance, bool verbose) {
00295 
00296   // Allocate vectors needed for iterations
00297   Epetra_Vector phat(x.Map()); phat.SetFlopCounter(x);
00298   Epetra_Vector p(x.Map()); p.SetFlopCounter(x);
00299   Epetra_Vector shat(x.Map()); shat.SetFlopCounter(x);
00300   Epetra_Vector s(x.Map()); s.SetFlopCounter(x);
00301   Epetra_Vector r(x.Map()); r.SetFlopCounter(x);
00302   Epetra_Vector rtilde(x.Map()); rtilde.Random(); rtilde.SetFlopCounter(x);
00303   Epetra_Vector v(x.Map()); v.SetFlopCounter(x);
00304 
00305 
00306   A.Multiply(false, x, r); // r = A*x
00307 
00308   r.Update(1.0, b, -1.0); // r = b - A*x
00309 
00310   double r_norm, b_norm, scaled_r_norm, rhon, rhonm1 = 1.0;
00311   double alpha = 1.0, omega = 1.0, sigma;
00312   double omega_num, omega_den;
00313   r.Norm2(&r_norm);
00314   b.Norm2(&b_norm);
00315   scaled_r_norm = r_norm/b_norm;
00316   r.Dot(rtilde,&rhon);
00317 
00318   if (verbose) cout << "Initial residual = " << r_norm
00319         << " Scaled residual = " << scaled_r_norm << endl;
00320 
00321 
00322   for (int i=0; i<Maxiter; i++) { // Main iteration loop
00323 
00324     double beta = (rhon/rhonm1) * (alpha/omega);
00325     rhonm1 = rhon;
00326 
00327     /* p    = r + beta*(p - omega*v)       */
00328     /* phat = M^-1 p                       */
00329     /* v    = A phat                       */
00330 
00331     double dtemp = - beta*omega;
00332 
00333     p.Update(1.0, r, dtemp, v, beta);
00334     if (M==0)
00335       phat.Scale(1.0, p);
00336     else
00337       M->Solve(false, p, phat);
00338     A.Multiply(false, phat, v);
00339 
00340 
00341     rtilde.Dot(v,&sigma);
00342     alpha = rhon/sigma;
00343 
00344     /* s = r - alpha*v                     */
00345     /* shat = M^-1 s                       */
00346     /* r = A shat (r is a tmp here for t ) */
00347 
00348     s.Update(-alpha, v, 1.0, r, 0.0);
00349     if (M==0)
00350       shat.Scale(1.0, s);
00351     else
00352       M->Solve(false, s, shat);
00353     A.Multiply(false, shat, r);
00354 
00355     r.Dot(s, &omega_num);
00356     r.Dot(r, &omega_den);
00357     omega = omega_num / omega_den;
00358 
00359     /* x = x + alpha*phat + omega*shat */
00360     /* r = s - omega*r */
00361 
00362     x.Update(alpha, phat, omega, shat, 1.0);
00363     r.Update(1.0, s, -omega);
00364 
00365     r.Norm2(&r_norm);
00366     scaled_r_norm = r_norm/b_norm;
00367     r.Dot(rtilde,&rhon);
00368 
00369     if (verbose) cout << "Iter "<< i << " residual = " << r_norm
00370           << " Scaled residual = " << scaled_r_norm << endl;
00371 
00372     if (scaled_r_norm < Tolerance) break;
00373   }
00374   return;
00375 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines